-
Dave Whipp authored
Added option for a pre-run isostatic compensation stage with a limit on the rate of compensation based on the geometry and a defined max velocity
Dave Whipp authoredAdded option for a pre-run isostatic compensation stage with a limit on the rate of compensation based on the geometry and a defined max velocity
define_bc.f90 26.26 KiB
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
! |
! ||===\\ |
! || \\ |
! || || //==\\ || || //==|| ||/==\\ |
! || || || || || || || || || || |
! || // || || || || || || || |
! ||===// \\==// \\==\\ \\==\\ || |
! |
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
! |
! DEFINE_BC Feb. 2009 |
! |
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
subroutine define_bc (params,osolve,vo,bcdef,nest)
!------------------------------------------------------------------------------|
!(((((((((((((((( Purpose of the routine )))))))))))))))))))))))))))))))))))))
!------------------------------------------------------------------------------|
! This defines the velocity boundary conditions
!------------------------------------------------------------------------------|
!(((((((((((((((( declaration of the subroutine arguments ))))))))))))))))))))
!------------------------------------------------------------------------------|
use definitions
implicit none
type (parameters) params
type (octreesolve) osolve
type (void) vo
type (bc_definition) bcdef
type (nest_info) nest
!------------------------------------------------------------------------------|
!(((((((((((((((( declaration of the subroutine internal variables )))))))))))))
!------------------------------------------------------------------------------|
integer i,iproc,nproc,ierr
double precision eps,lsf0,pi,lorig,h,x1,x2,phi,spu,spv,yend,cper,cscl,xstart
double precision ystart,theta,l,vin,vzfluxscl,cntvel,dxy,xend,xsym,ymax,xwidth
double precision ywidth
double precision,dimension(:),allocatable :: x0,ldisp
integer ie,ij,j,jp,nelemx,nelemz
logical firstcall
character(len=40) :: bccase
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
INCLUDE 'mpif.h'
call mpi_comm_size (mpi_comm_world,nproc,ierr)
call mpi_comm_rank (mpi_comm_world,iproc,ierr)
eps=1.d-10
osolve%kfix=0
osolve%kfixt=0
if (bcdef%bctype=='') then
bccase=params%infile
else
bccase=bcdef%bctype
endif
if (params%nest) bccase='nest'
select case(trim(bccase))
case('input.txt','input.small.txt','basic')
if (trim(bcdef%bcorder).ne.'xyz' .and. trim(bcdef%bcorder).ne.'xzy' .and. &
trim(bcdef%bcorder).ne.'yxz' .and. trim(bcdef%bcorder).ne.'yzx' .and. &
trim(bcdef%bcorder).ne.'zxy' .and. trim(bcdef%bcorder).ne.'zyx') then
if (iproc==0) write (*,*) 'No/unsupported BC order provided, assuming bcorder=xyz'
endif
do i=1,osolve%nnode
select case(trim(bcdef%bcorder))
case ('xyz')
if (osolve%x(i).lt.eps) then
if (bcdef%fixux0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux0
endif
if (bcdef%fixvx0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx0
endif
if (bcdef%fixwx0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx0
endif
endif
if (osolve%x(i).gt.1.d0-eps) then
if (bcdef%fixux1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux1
endif
if (bcdef%fixvx1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx1
endif
if (bcdef%fixwx1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx1
endif
endif
if (osolve%y(i).lt.eps) then
if (bcdef%fixuy0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy0
endif
if (bcdef%fixvy0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy0
endif
if (bcdef%fixwy0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy0
endif
endif
if (osolve%y(i).gt.1.d0-eps) then
if (bcdef%fixuy1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy1
endif
if (bcdef%fixvy1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy1
endif
if (bcdef%fixwy1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy1
endif
endif
if (osolve%z(i).lt.eps) then
if (bcdef%fixuz0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz0
endif
if (bcdef%fixvz0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz0
endif
if (bcdef%fixwz0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz0
endif
osolve%kfixt(i)=1
osolve%temp(i)=1.d0
endif
if (osolve%z(i).gt.1.d0-eps) then
if (bcdef%fixuz1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz1
endif
if (bcdef%fixvz1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz1
endif
if (bcdef%fixwz1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz1
endif
osolve%kfixt(i)=1
osolve%temp(i)=0.d0
endif
case('xzy')
if (osolve%x(i).lt.eps) then
if (bcdef%fixux0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux0
endif
if (bcdef%fixvx0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx0
endif
if (bcdef%fixwx0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx0
endif
endif
if (osolve%x(i).gt.1.d0-eps) then
if (bcdef%fixux1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux1
endif
if (bcdef%fixvx1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx1
endif
if (bcdef%fixwx1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx1
endif
endif
if (osolve%z(i).lt.eps) then
if (bcdef%fixuz0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz0
endif
if (bcdef%fixvz0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz0
endif
if (bcdef%fixwz0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz0
endif
osolve%kfixt(i)=1
osolve%temp(i)=1.d0
endif
if (osolve%z(i).gt.1.d0-eps) then
if (bcdef%fixuz1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz1
endif
if (bcdef%fixvz1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz1
endif
if (bcdef%fixwz1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz1
endif
osolve%kfixt(i)=1
osolve%temp(i)=0.d0
endif
if (osolve%y(i).lt.eps) then
if (bcdef%fixuy0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy0
endif
if (bcdef%fixvy0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy0
endif
if (bcdef%fixwy0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy0
endif
endif
if (osolve%y(i).gt.1.d0-eps) then
if (bcdef%fixuy1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy1
endif
if (bcdef%fixvy1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy1
endif
if (bcdef%fixwy1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy1
endif
endif
case('yxz')
if (osolve%y(i).lt.eps) then
if (bcdef%fixuy0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy0
endif
if (bcdef%fixvy0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy0
endif
if (bcdef%fixwy0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy0
endif
endif
if (osolve%y(i).gt.1.d0-eps) then
if (bcdef%fixuy1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy1
endif
if (bcdef%fixvy1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy1
endif
if (bcdef%fixwy1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy1
endif
endif
if (osolve%x(i).lt.eps) then
if (bcdef%fixux0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux0
endif
if (bcdef%fixvx0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx0
endif
if (bcdef%fixwx0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx0
endif
endif
if (osolve%x(i).gt.1.d0-eps) then
if (bcdef%fixux1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux1
endif
if (bcdef%fixvx1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx1
endif
if (bcdef%fixwx1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx1
endif
endif
if (osolve%z(i).lt.eps) then
if (bcdef%fixuz0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz0
endif
if (bcdef%fixvz0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz0
endif
if (bcdef%fixwz0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz0
endif
osolve%kfixt(i)=1
osolve%temp(i)=1.d0
endif
if (osolve%z(i).gt.1.d0-eps) then
if (bcdef%fixuz1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz1
endif
if (bcdef%fixvz1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz1
endif
if (bcdef%fixwz1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz1
endif
osolve%kfixt(i)=1
osolve%temp(i)=0.d0
endif
case('yzx')
if (osolve%y(i).lt.eps) then
if (bcdef%fixuy0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy0
endif
if (bcdef%fixvy0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy0
endif
if (bcdef%fixwy0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy0
endif
endif
if (osolve%y(i).gt.1.d0-eps) then
if (bcdef%fixuy1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy1
endif
if (bcdef%fixvy1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy1
endif
if (bcdef%fixwy1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy1
endif
endif
if (osolve%z(i).lt.eps) then
if (bcdef%fixuz0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz0
endif
if (bcdef%fixvz0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz0
endif
if (bcdef%fixwz0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz0
endif
osolve%kfixt(i)=1
osolve%temp(i)=1.d0
endif
if (osolve%z(i).gt.1.d0-eps) then
if (bcdef%fixuz1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz1
endif
if (bcdef%fixvz1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz1
endif
if (bcdef%fixwz1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz1
endif
osolve%kfixt(i)=1
osolve%temp(i)=0.d0
endif
if (osolve%x(i).lt.eps) then
if (bcdef%fixux0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux0
endif
if (bcdef%fixvx0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx0
endif
if (bcdef%fixwx0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx0
endif
endif
if (osolve%x(i).gt.1.d0-eps) then
if (bcdef%fixux1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux1
endif
if (bcdef%fixvx1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx1
endif
if (bcdef%fixwx1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx1
endif
endif
case('zxy')
if (osolve%z(i).lt.eps) then
if (bcdef%fixuz0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz0
endif
if (bcdef%fixvz0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz0
endif
if (bcdef%fixwz0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz0
endif
osolve%kfixt(i)=1
osolve%temp(i)=1.d0
endif
if (osolve%z(i).gt.1.d0-eps) then
if (bcdef%fixuz1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz1
endif
if (bcdef%fixvz1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz1
endif
if (bcdef%fixwz1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz1
endif
osolve%kfixt(i)=1
osolve%temp(i)=0.d0
endif
if (osolve%x(i).lt.eps) then
if (bcdef%fixux0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux0
endif
if (bcdef%fixvx0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx0
endif
if (bcdef%fixwx0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx0
endif
endif
if (osolve%x(i).gt.1.d0-eps) then
if (bcdef%fixux1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux1
endif
if (bcdef%fixvx1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx1
endif
if (bcdef%fixwx1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx1
endif
endif
if (osolve%y(i).lt.eps) then
if (bcdef%fixuy0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy0
endif
if (bcdef%fixvy0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy0
endif
if (bcdef%fixwy0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy0
endif
endif
if (osolve%y(i).gt.1.d0-eps) then
if (bcdef%fixuy1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy1
endif
if (bcdef%fixvy1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy1
endif
if (bcdef%fixwy1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy1
endif
endif
case('zyx')
if (osolve%z(i).lt.eps) then
if (bcdef%fixuz0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz0
endif
if (bcdef%fixvz0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz0
endif
if (bcdef%fixwz0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz0
endif
osolve%kfixt(i)=1
osolve%temp(i)=1.d0
endif
if (osolve%z(i).gt.1.d0-eps) then
if (bcdef%fixuz1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz1
endif
if (bcdef%fixvz1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz1
endif
if (bcdef%fixwz1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz1
endif
osolve%kfixt(i)=1
osolve%temp(i)=0.d0
endif
if (osolve%y(i).lt.eps) then
if (bcdef%fixuy0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy0
endif
if (bcdef%fixvy0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy0
endif
if (bcdef%fixwy0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy0
endif
endif
if (osolve%y(i).gt.1.d0-eps) then
if (bcdef%fixuy1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy1
endif
if (bcdef%fixvy1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy1
endif
if (bcdef%fixwy1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy1
endif
endif
if (osolve%x(i).lt.eps) then
if (bcdef%fixux0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux0
endif
if (bcdef%fixvx0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx0
endif
if (bcdef%fixwx0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx0
endif
endif
if (osolve%x(i).gt.1.d0-eps) then
if (bcdef%fixux1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux1
endif
if (bcdef%fixvx1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx1
endif
if (bcdef%fixwx1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx1
endif
endif
case default
if (osolve%x(i).lt.eps) then
if (bcdef%fixux0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux0
endif
if (bcdef%fixvx0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx0
endif
if (bcdef%fixwx0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx0
endif
endif
if (osolve%x(i).gt.1.d0-eps) then
if (bcdef%fixux1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%ux1
endif
if (bcdef%fixvx1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vx1
endif
if (bcdef%fixwx1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wx1
endif
endif
if (osolve%y(i).lt.eps) then
if (bcdef%fixuy0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy0
endif
if (bcdef%fixvy0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy0
endif
if (bcdef%fixwy0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy0
endif
endif
if (osolve%y(i).gt.1.d0-eps) then
if (bcdef%fixuy1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uy1
endif
if (bcdef%fixvy1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vy1
endif
if (bcdef%fixwy1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wy1
endif
endif
if (osolve%z(i).lt.eps) then
if (bcdef%fixuz0) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz0
endif
if (bcdef%fixvz0) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz0
endif
if (bcdef%fixwz0) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz0
endif
osolve%kfixt(i)=1
osolve%temp(i)=1.d0
endif
if (osolve%z(i).gt.1.d0-eps) then
if (bcdef%fixuz1) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=bcdef%uz1
endif
if (bcdef%fixvz1) then
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=bcdef%vz1
endif
if (bcdef%fixwz1) then
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=bcdef%wz1
endif
osolve%kfixt(i)=1
osolve%temp(i)=0.d0
endif
end select
if (.not.vo%influid(i)) then
osolve%kfixt(i)=1
osolve%temp(i)=0.d0
endif
if (abs(bcdef%utrans).gt.eps) then
osolve%u(i)=osolve%u(i)+bcdef%utrans
endif
if (abs(bcdef%vtrans).gt.eps) then
osolve%v(i)=osolve%v(i)+bcdef%vtrans
endif
enddo
case('iso_only')
do i=1,osolve%nnode
if (osolve%x(i).lt.eps) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=0.d0
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=0.d0
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=0.d0
endif
if (osolve%x(i).gt.1.d0-eps) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=0.d0
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=0.d0
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=0.d0
endif
if (osolve%y(i).lt.eps) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=0.d0
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=0.d0
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=0.d0
endif
if (osolve%y(i).gt.1.d0-eps) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=0.d0
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=0.d0
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=0.d0
endif
if (osolve%z(i).lt.eps) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=0.d0
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=0.d0
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=0.d0
osolve%kfixt(i)=1
osolve%temp(i)=1.d0
endif
if (osolve%z(i).gt.1.d0-eps) then
osolve%kfix((i-1)*3+1)=1
osolve%u(i)=0.d0
osolve%kfix((i-1)*3+2)=1
osolve%v(i)=0.d0
osolve%kfix((i-1)*3+3)=1
osolve%w(i)=0.d0
osolve%kfixt(i)=1
osolve%temp(i)=0.d0
endif
if (.not.vo%influid(i)) then
osolve%kfixt(i)=1
osolve%temp(i)=0.d0
endif
enddo
case ('input.2Dpunch','2Dpunch')
call define_bc_2Dpunch (osolve%nnode,osolve%kfix,osolve%kfixt,osolve%x,osolve%y,osolve%z,osolve%u,osolve%v,osolve%w,osolve%temp,vo)
case ('input.2Dpunch_vert','2Dpunch_vert')
call define_bc_2Dpunch_vert (osolve%nnode,osolve%kfix,osolve%kfixt,osolve%x,osolve%y,osolve%z,osolve%u,osolve%v,osolve%w,osolve%temp,vo)
case ('input.3Dpunch','3Dpunch')
call define_bc_3Dpunch (osolve%nnode,osolve%kfix,osolve%kfixt,osolve%x,osolve%y,osolve%z,osolve%u,osolve%v,osolve%w,osolve%temp,vo)
case ('input.folding','folding')
call define_bc_folding (osolve%nnode,osolve%kfix,osolve%kfixt,osolve%x,osolve%y,osolve%z,osolve%u,osolve%v,osolve%w,osolve%temp,vo)
case ('input.jgr','jgr')
call define_bc_jgr (osolve%nnode,osolve%kfix,osolve%kfixt,osolve%x,osolve%y,osolve%z,osolve%u,osolve%v,osolve%w,osolve%temp,vo)
case ('input.nest','nest')
call define_bc_nest (osolve%nnode,osolve%kfix,osolve%kfixt,osolve%x,osolve%y,osolve%z,osolve%u,osolve%v,osolve%w,osolve%temp,vo,params,nest)
case ('input.parallipipede','parallipipede')
call define_bc_parallipipede (osolve%nnode,osolve%kfix,osolve%kfixt,osolve%x,osolve%y,osolve%z,osolve%u,osolve%v,osolve%w,osolve%temp,vo)
case ('input.pipo','pipo')
call define_bc_pipo (osolve%nnode,osolve%kfix,osolve%kfixt,osolve%x,osolve%y,osolve%z,osolve%u,osolve%v,osolve%w,osolve%temp,vo)
case ('input.riedel','riedel')
call define_bc_riedel (osolve%nnode,osolve%kfix,osolve%kfixt,osolve%x,osolve%y,osolve%z,osolve%u,osolve%v,osolve%w,osolve%temp,vo)
case ('input.ritske','ritske')
call define_bc_ritske (osolve%nnode,osolve%kfix,osolve%kfixt,osolve%x,osolve%y,osolve%z,osolve%u,osolve%v,osolve%w,osolve%temp,vo)
case ('input.ritske_isurf','ritske_isurf')
call define_bc_ritske_isurf (osolve%nnode,osolve%kfix,osolve%kfixt,osolve%x,osolve%y,osolve%z,osolve%u,osolve%v,osolve%w,osolve%temp,vo)
case ('input.segmented_s_line','segmented_s_line')
call define_bc_segmented_s_line (params,osolve,vo,bcdef,nest)
case ('input.sphere','sphere')
call define_bc_sphere (osolve%nnode,osolve%kfix,osolve%kfixt,osolve%x,osolve%y,osolve%z,osolve%u,osolve%v,osolve%w,osolve%temp,vo)
case ('input.subduction','subduction')
call define_bc_subduction (osolve%nnode,osolve%kfix,osolve%kfixt,osolve%x,osolve%y,osolve%z,osolve%u,osolve%v,osolve%w,osolve%temp,vo)
case default
if (iproc.eq.0) print *,params%infile
call stop_run (' pb with input file in define_bc.f90')
end select
! here we need to impose that the bc satisfy the bad faces constrains
do ie=1,osolve%nface
do j=1,4
jp=1+mod(j,4)
i=osolve%iface(j+4,ie)
if (osolve%kfix((i-1)*3+1).eq.1) osolve%u(i)=(osolve%u(osolve%iface(j,ie))+osolve%u(osolve%iface(jp,ie)))/2.d0
if (osolve%kfix((i-1)*3+2).eq.1) osolve%v(i)=(osolve%v(osolve%iface(j,ie))+osolve%v(osolve%iface(jp,ie)))/2.d0
if (osolve%kfix((i-1)*3+3).eq.1) osolve%w(i)=(osolve%w(osolve%iface(j,ie))+osolve%w(osolve%iface(jp,ie)))/2.d0
if (osolve%kfixt(i).eq.1) osolve%temp(i)=(osolve%temp(osolve%iface(j,ie))+osolve%temp(osolve%iface(jp,ie)))/2.d0
enddo
i=osolve%iface(9,ie)
if (osolve%kfix((i-1)*3+1).eq.1) osolve%u(i)=sum(osolve%u(osolve%iface(1:4,ie)))/4.d0
if (osolve%kfix((i-1)*3+2).eq.1) osolve%v(i)=sum(osolve%v(osolve%iface(1:4,ie)))/4.d0
if (osolve%kfix((i-1)*3+3).eq.1) osolve%w(i)=sum(osolve%w(osolve%iface(1:4,ie)))/4.d0
if (osolve%kfixt(i).eq.1) osolve%temp(i)=sum(osolve%temp(osolve%iface(1:4,ie)))/4.d0
enddo
if (params%debug==2) call output_bc (osolve)
return
end
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|