!
! This file is part of SACAMOS, State of the Art CAble MOdels in Spice.
! It was developed by the University of Nottingham and the Netherlands Aerospace
! Centre (NLR) for ESA under contract number 4000112765/14/NL/HK.
!
! Copyright (C) 2016-2017 University of Nottingham
!
! SACAMOS is free software: you can redistribute it and/or modify it under the
! terms of the GNU General Public License as published by the Free Software
! Foundation, either version 3 of the License, or (at your option) any later
! version.
!
! SACAMOS is distributed in the hope that it will be useful, but
! WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
! or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
! for more details.
!
! A copy of the GNU General Public License version 3 can be found in the
! file GNU_GPL_v3 in the root or at .
!
! SACAMOS uses the EISPACK library (in /SRC/EISPACK). EISPACK is subject to
! the GNU Lesser General Public License. A copy of the GNU Lesser General Public
! License version can be found in the file GNU_LGPL in the root of EISPACK
! (/SRC/EISPACK ) or at .
!
! The University of Nottingham can be contacted at: ggiemr@nottingham.ac.uk
!
! File Contents:
! SUBROUTINE Dconnector_set_parameters
! SUBROUTINE Dconnector_set_internal_domain_information
! SUBROUTINE Dconnector_plot
!
! NAME
! Dconnector_set_parameters
!
! AUTHORS
! Chris Smartt
!
! DESCRIPTION
! Set the overall parameters for a Dconnector cable
!
! COMMENTS
!
!
! HISTORY
!
! started 12/4/2016 CJS
! 20/5/2016 CJS Allow use of the Laplace solver in the internaln domain
!
!
SUBROUTINE Dconnector_set_parameters(cable)
USE type_specifications
IMPLICIT NONE
! variables passed to subroutine
type(cable_specification_type),intent(INOUT) :: cable
! local variables
! START
cable%cable_type=cable_geometry_type_Dconnector
cable%tot_n_conductors=0 ! this is set in the cable specification file
cable%tot_n_domains=2
cable%n_external_conductors=1
cable%n_internal_conductors=0 ! this is set in the cable specification file
cable%n_internal_domains=1
cable%n_parameters=4
cable%n_dielectric_filters=0
cable%n_transfer_impedance_models=0
END SUBROUTINE Dconnector_set_parameters
!
! NAME
! Dconnector_set_internal_domain_information
!
! AUTHORS
! Chris Smartt
!
! DESCRIPTION
! Set the overall parameters for a Dconnector cable
!
! COMMENTS
! Set the dimension of the domain transformation matrices to include an external reference conductor for the cable
!
!
! HISTORY
!
! started 12/4/201 CJS
! 8/5/2017 CJS: Include references to Theory_Manual
!
!
SUBROUTINE Dconnector_set_internal_domain_information(cable)
USE type_specifications
USE constants
USE general_module
USE maths
USE PUL_parameter_module
IMPLICIT NONE
! variables passed to subroutine
type(cable_specification_type),intent(INOUT) :: cable
! local variables
integer :: nc
integer :: dim
integer :: domain
type(PUL_type) :: PUL
integer :: ierr
! variables for cable parameter checks
logical :: cable_spec_error
real(dp) :: rw
real(dp) :: p
real(dp) :: s
real(dp) :: o
! variables for cross section geometry
integer :: ncrow(2) ! number of conductors in the two rows
real(dp) :: w(2) ! width of the two rows of conductors
integer :: i
character(LEN=error_message_length) :: message
! START
! Check the cable parameters
nc=cable%tot_n_conductors
rw=cable%parameters(1)
p=cable%parameters(2)
s=cable%parameters(3)
o=cable%parameters(4)
cable_spec_error=.FALSE. ! assume no errors initially
message=''
! we must have 5 or more conductors in a D connector
if (nc.LT.5) then
message='ERROR: There must be 5 or more conductors in a D connector'
write(*,'(A)')message
cable_spec_error=.TRUE.
end if
CALL cylindrical_check(rw,cable_spec_error,cable%cable_name,message)
CALL gt_zero_check(p,cable_spec_error,cable%cable_name,message)
CALL gt_zero_check(s,cable_spec_error,cable%cable_name,message)
CALL gt_zero_check(o,cable_spec_error,cable%cable_name,message)
if (cable_spec_error) then
run_status='ERROR in cable_model_builder, error on parameters for cable:'//trim(cable%cable_name)//'. '//trim(message)
CALL write_program_status()
STOP 1
end if
! set the information related to the number of conductors.
! Note that the number of conductors includes the shield which encloses the connector pins
cable%n_internal_conductors=nc-1
cable%n_external_conductors=1
! Set the parameters for the single internal domain
domain=1
cable%n_internal_conductors_in_domain(domain)=nc
! The number of modes in the internal domain is nc-1
dim=nc-1
cable%L_domain(domain)%dim=dim
ALLOCATE(cable%L_domain(domain)%mat(dim,dim))
cable%C_domain(domain)%dim=dim
ALLOCATE(cable%C_domain(domain)%mat(dim,dim))
if (use_laplace) then
! allocate memory for the PUL parameter solver interface
if(verbose) write(*,*)'Domain:',domain
if(verbose) write(*,*)'Allocating PUL data structure for Dconnector'
CALL allocate_and_reset_PUL_data(PUL,nc)
PUL%shape(1:nc-1)=circle
! work out the number of conductors in the two rows
ncrow(1)=nc/2
ncrow(2)=(nc-1)-ncrow(1)
write(*,*)'Connector model: nc1=',ncrow(1),' nc2=',ncrow(2)
! width of the two rows of conductors
w(1)=(ncrow(1)-1)*p
w(2)=(ncrow(2)-1)*p
! top row conductors
do i=1,ncrow(1)
PUL%x(i)=-w(1)/2d0+(i-1)*p
PUL%y(i)=s/2d0
PUL%r(i)=rw
PUL%rd(i)=rw
PUL%epsr(i)=1d0 ! no dielectric on conductors, they are in a homogeneous dielectric medium
end do
! bottom row conductors
do i=1,ncrow(2)
PUL%x(ncrow(1)+i)=-w(2)/2d0+(i-1)*p
PUL%y(ncrow(1)+i)=-s/2d0
PUL%r(ncrow(1)+i)=rw
PUL%rd(ncrow(1)+i)=rw
PUL%epsr(ncrow(1)+i)=1d0 ! no dielectric on conductors, they are in a homogeneous dielectric medium
end do
PUL%epsr_background = 1d0 ! permittivity of homogeneous dielectric medium surrounding conductors
! no ground plane
PUL%ground_plane_present=.FALSE.
! add overshield i.e. the Dconnector shield
PUL%overshield_present=.TRUE.
PUL%overshield_shape=Dshape
PUL%overshield_x = 0d0 ! shield is centred at the origin in this calculation
PUL%overshield_y = 0d0
PUL%overshield_r = rw+o ! Dconnector curve radius (dl is defined in terms of this)
PUL%overshield_w = w(1)+2d0*(rw+o) ! Dconnector shield top width
PUL%overshield_w2 = w(2)+2d0*(rw+o) ! Dconnector shield bottom width
PUL%overshield_h = s+2d0*(rw+o) ! Dconnector shield height
CALL PUL_LC_Laplace(PUL,cable%cable_name,cable%Y_fit_model_order,cable%Y_fit_freq_spec,domain)
cable%L_domain(domain)%mat(:,:)=PUL%L%mat(:,:)
cable%C_domain(domain)%mat(:,:)=PUL%C%mat(:,:)
else
run_status='ERROR in Dconnector: this cable type must use the Laplace solver for the PUL paramter calculation'
CALL write_program_status()
STOP 1
end if
CALL Z_Y_from_L_C(cable%L_domain(domain),cable%C_domain(domain),cable%Z_domain(domain),cable%Y_domain(domain))
if (use_laplace) CALL deallocate_PUL_data(PUL) ! deallocate the PUL data structure
! Set the domain decomposition matrices ! Theory_Manual_Eqn 6.19, 6.20
! The dimension of the domain transformation matrices is the number of conductors+1
dim=nc+1
cable%MI%dim=dim
ALLOCATE(cable%MI%mat(dim,dim))
cable%MI%mat(:,:)=0d0
cable%MV%dim=dim
ALLOCATE(cable%MV%mat(dim,dim))
cable%MV%mat(:,:)=0d0
do i=1,nc-1 ! loop over conductors except the local reference (shield) conductor
cable%MI%mat(i,i)=1d0
cable%MV%mat(i,i)= 1d0
cable%MV%mat(i,nc)=-1d0
end do
! row for local reference conductor
do i=1,nc
cable%MI%mat(nc,i)=1d0
end do
cable%MV%mat(nc,nc)=1d0
cable%MV%mat(nc,nc+1)=-1d0
! row for global reference
do i=1,nc+1
cable%MI%mat(nc+1,i)=1d0
end do
cable%MV%mat(nc+1,nc+1)=1d0
! Set the local reference conductor numbering
ALLOCATE( cable%local_reference_conductor(nc) )
do i=1,nc-1
cable%local_reference_conductor(i)=nc ! inner conductor, reference is the shield conductor
end do
cable%local_reference_conductor(nc)=0 ! external domain conductor, reference not known
! Set the local domain information: include a reference conductor in the count
ALLOCATE( cable%local_domain_n_conductors(1:cable%tot_n_domains) )
cable%local_domain_n_conductors(1)=nc ! internal domain
cable%local_domain_n_conductors(2)=2 ! external domain
! Set the external domain conductor and dielectric information
ALLOCATE( cable%external_model(cable%n_external_conductors) )
CALL reset_external_conductor_model(cable%external_model(1))
cable%external_model(1)%conductor_type=Dshape
cable%external_model(1)%conductor_radius=o+rw
cable%external_model(1)%conductor_width=w(1)+2d0*(o+rw)
cable%external_model(1)%conductor_width2=w(2)+2d0*(o+rw)
cable%external_model(1)%conductor_height=s+2d0*(o+rw)
END SUBROUTINE Dconnector_set_internal_domain_information
!
! NAME
! Dconnector_plot
!
! AUTHORS
! Chris Smartt
!
! DESCRIPTION
! plot Dconnector cable
!
! COMMENTS
! the angle has an impact here...
! The conductor geometry must be consistent with the documentation...
!
! HISTORY
!
! started 14/4/2016 CJS
!
!
SUBROUTINE Dconnector_plot(cable,x_offset,y_offset,theta,xmin,xmax,ymin,ymax)
USE type_specifications
USE general_module
IMPLICIT NONE
! variables passed to subroutine
type(cable_specification_type),intent(IN) :: cable
real(dp),intent(IN) :: x_offset,y_offset,theta
real(dp),intent(INOUT) :: xmin,xmax,ymin,ymax
! local variables
! variables for cable parameters
integer :: nc
real(dp) :: rw
real(dp) :: p
real(dp) :: s
real(dp) :: o
! variables for cross section geometry
integer :: ncrow(2) ! number of conductors in the two rows
real(dp) :: w(2) ! width of the two rows of conductors
real(dp) :: cx,cy,x,y
real(dp) :: wd,wd2,hd
integer :: i
! START
nc=cable%tot_n_conductors
rw=cable%parameters(1)
p=cable%parameters(2)
s=cable%parameters(3)
o=cable%parameters(4)
! plot inner conductors
! work out the number of conductors in the two rows
ncrow(1)=nc/2
ncrow(2)=(nc-1)-ncrow(1)
! width of the two rows of conductors
w(1)=(ncrow(1)-1)*p
w(2)=(ncrow(2)-1)*p
! top row conductors
do i=1,ncrow(1)
cx=-w(1)/2d0+(i-1)*p
cy=s/2d0
x=x_offset+cx*cos(theta)-cy*sin(theta)
y=y_offset+cx*sin(theta)+cy*cos(theta)
CALL write_circle(x,y,rw,conductor_geometry_file_unit,xmin,xmax,ymin,ymax)
end do
! bottom row conductors
do i=1,ncrow(2)
cx=-w(2)/2d0+(i-1)*p
cy=-s/2d0
x=x_offset+cx*cos(theta)-cy*sin(theta)
y=y_offset+cx*sin(theta)+cy*cos(theta)
CALL write_circle(x,y,rw,conductor_geometry_file_unit,xmin,xmax,ymin,ymax)
end do
! write the outer conductor
x=x_offset
y=y_offset
CALL write_Dshape(x,y,w(1)/2d0,w(2)/2d0,s/2d0,rw+o,theta,conductor_geometry_file_unit,xmin,xmax,ymin,ymax)
RETURN
END SUBROUTINE Dconnector_plot