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_tel3(ielfa,xrlfa,icyclic,ifatie,gradtfa, 20 & gradtel,c,ipnei,xxi,ifabou,xxn,xload,nfacea,nfaceb,ncfd) 21! 22! interpolate/extrapolate the temperature gradient from the 23! center of the elements to the center of the faces 24! 25 implicit none 26! 27 integer ielfa(4,*),icyclic,ifatie(*),ipnei(*),nfacea,nfaceb, 28 & iel1,iel2,i,l,indexf,ifabou(*),ncfd 29! 30 real*8 xrlfa(3,*),gradtfa(3,*),gradtel(3,*),c(3,3),xxi(3,*), 31 & gradnor,xl1,xl2,q,xxn(3,*),xload(2,*) 32! 33! 34! 35 do i=nfacea,nfaceb 36 iel1=ielfa(1,i) 37 xl1=xrlfa(1,i) 38 iel2=ielfa(2,i) 39 if(iel2.gt.0) then 40! 41! face between two elements 42! 43 xl2=xrlfa(2,i) 44 if((icyclic.eq.0).or.(ifatie(i).eq.0)) then 45 do l=1,ncfd 46 gradtfa(l,i)=xl1*gradtel(l,iel1)+ 47 & xl2*gradtel(l,iel2) 48 enddo 49 elseif(ifatie(i).gt.0) then 50 do l=1,ncfd 51 gradtfa(l,i)=xl1*gradtel(l,iel1)+xl2* 52 & (gradtel(1,iel2)*c(l,1)+ 53 & gradtel(2,iel2)*c(l,2)+ 54 & gradtel(3,iel2)*c(l,3)) 55 enddo 56 else 57 do l=1,ncfd 58 gradtfa(l,i)=xl1*gradtel(l,iel1)+xl2* 59 & (gradtel(1,iel2)*c(1,l)+ 60 & gradtel(2,iel2)*c(2,l)+ 61 & gradtel(3,iel2)*c(3,l)) 62 enddo 63 endif 64 elseif(ielfa(3,i).gt.0) then 65! 66! the above condition implies iel2!=0 67! if iel2 were zero, no b.c. would apply and 68! ielfa(3,i) would be zero or negative 69! 70! boundary face; no zero gradient 71! 72 if(ifabou(-iel2).gt.0) then 73! 74! temperature given: extrapolate gradient 75! 76 do l=1,ncfd 77 gradtfa(l,i)=xl1*gradtel(l,iel1)+ 78 & xrlfa(3,i)*gradtel(l,abs(ielfa(3,i))) 79 enddo 80 else 81! 82! facial temperature gradient given, may be nonzero 83! 84 indexf=ipnei(iel1)+ielfa(4,i) 85! 86 if(ifabou(-iel2+6).eq.0) then 87 q=0.d0 88 else 89 q=xload(1,ifabou(-iel2+6)) 90 endif 91 gradnor=gradtel(1,iel1)*xxn(1,indexf) 92 & +gradtel(2,iel1)*xxn(2,indexf) 93 & +gradtel(3,iel1)*xxn(3,indexf)-q 94 do l=1,ncfd 95 gradtfa(l,i)=gradtel(l,iel1) 96 & -gradnor*xxn(l,indexf) 97 enddo 98 endif 99 else 100! 101! boundary face; zero gradient 102! 103 indexf=ipnei(iel1)+ielfa(4,i) 104 gradnor=gradtel(1,iel1)*xxi(1,indexf)+ 105 & gradtel(2,iel1)*xxi(2,indexf)+ 106 & gradtel(3,iel1)*xxi(3,indexf) 107 do l=1,ncfd 108 gradtfa(l,i)=gradtel(l,iel1) 109 & -gradnor*xxi(l,indexf) 110 enddo 111 endif 112 enddo 113! 114 return 115 end 116