1! 2! CalculiX - A 3-dimensional finite element program 3! Copyright (C) 1998-2021 Guido Dhondt 4! 5! This program is free software; you can redistribute it and/or 6! modify it under the terms of the GNU General Public License as 7! published by the Free Software Foundation(version 2); 8! 9! 10! This program is distributed in the hope that it will be useful, 11! but WITHOUT ANY WARRANTY; without even the implied warranty of 12! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 13! GNU General Public License for more details. 14! 15! You should have received a copy of the GNU General Public License 16! along with this program; if not, write to the Free Software 17! Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA. 18! 19 subroutine extrapol_oel1(ielfa,xrlfa,vfap,vel,ifabou, 20 & nef,umfa,constant,dy,vfa,inlet,nfacea,nfaceb) 21! 22! extrapolation of turbulent kinetic energy values to the faces 23! 24 implicit none 25! 26 integer ielfa(4,*),ifabou(*),nfacea,nfaceb,nef,i,iel1,iel2, 27 & ipointer,inlet(*) 28! 29 real*8 xrlfa(3,*),vfap(0:7,*),vel(nef,0:7),xl1,dy(*), 30 & umfa(*),constant,vfa(0:7,*) 31! 32! 33! 34 do i=nfacea,nfaceb 35 iel1=ielfa(1,i) 36 xl1=xrlfa(1,i) 37 iel2=ielfa(2,i) 38 if(iel2.gt.0) then 39! 40! face between two elements: interpolation 41! 42 vfap(7,i)=xl1*vel(iel1,7)+xrlfa(2,i)*vel(iel2,7) 43! 44 elseif(ielfa(3,i).gt.0) then 45! 46! boundary face; no zero gradient 47! 48! iel2=0 is not possible: if iel2=0, there are no 49! boundary conditions on the face, hence it is an 50! exit, which means zero gradient and 51! ielfa(3,i) <= 0 52! 53 ipointer=-iel2 54! 55 if(ifabou(ipointer+5).gt.0) then 56! 57! wall: turbulent dissipation rate known 58! 59 vfap(7,i)=dy(ifabou(ipointer+5))*umfa(i)/vfa(5,i) 60 elseif((inlet(i).eq.1).or. 61 & (ifabou(ipointer+5).lt.0)) then 62! 63! inlet or sliding conditions: kinetic turbulent energy known 64! 65 vfap(7,i)=constant*umfa(i) 66 else 67! 68! extrapolation 69! 70 vfap(7,i)=xl1*vel(iel1,7)+xrlfa(3,i)*vel(ielfa(3,i),7) 71 endif 72 else 73! 74! boundary face; zero gradient in i-direction 75! 76 vfap(7,i)=vel(iel1,7) 77 endif 78 enddo 79! 80 return 81 end 82