!
! 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 shielded_twisted_pair_set_parameters
! SUBROUTINE shielded_twisted_pair_set_internal_domain_information
! SUBROUTINE shielded_twisted_pair_plot!
! NAME
! shielded_twisted_pair_set_parameters
!
! AUTHORS
! Chris Smartt
!
! DESCRIPTION
! Set the overall parameters for a shielded_twisted_pair cable
!
! COMMENTS
!
!
! HISTORY
!
! started 5/9/2016 CJS based on shielded_twisted_pair.F90
! 16/11/2017 CJS Include network synthesis process to replace s-domain transfer functions
!
!
SUBROUTINE shielded_twisted_pair_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_shielded_twisted_pair
cable%tot_n_conductors=3
cable%tot_n_domains=3
cable%n_external_conductors=1
cable%n_internal_conductors=2
cable%n_internal_domains=2 ! note we have both a common mode and a differential mode internal domain
cable%n_parameters=8
cable%n_dielectric_filters=2
cable%n_transfer_impedance_models=1
END SUBROUTINE shielded_twisted_pair_set_parameters
!
! NAME
! shielded_twisted_pair_set_internal_domain_information
!
! AUTHORS
! Chris Smartt
!
! DESCRIPTION
! Set the overall parameters for a shielded_twisted_pair cable
!
! COMMENTS
! Set the dimension of the domain transformation matrices to include an external reference conductor for the cable
!
!
! HISTORY
!
! started 5/9/2016 CJS based on shielded_twisted_pair.F90
! 8/9/2016 CJS common mode/ differential mode loss correction
! 19/9/2016 CJS frequency dependent dielectric in Laplace solver
! 2/11/2016 CJS inhomogeneous dielectric in twisted pair model
! 15/12/2016 CJS frequeny dependent dielectric in twisted pair model with Laplace solver
! 8/5/2017 CJS: Include references to Theory_Manual
!
!
SUBROUTINE shielded_twisted_pair_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 :: n_conductors
integer :: dim
integer :: domain
type(PUL_type) :: PUL
real(dp) :: L11,L12
real(dp) :: C11,C12
real(dp) :: LC,LD,CC,CD
logical :: dielectric_is_homogeneous
real(dp) :: C_air
type(Sfilter) :: jw
integer :: i
real(dp) :: epsr
! variables for cable parameter checks
logical :: cable_spec_error
real(dp) :: rw
real(dp) :: s
real(dp) :: rd
real(dp) :: rs
real(dp) :: rd2
real(dp) :: t
real(dp) :: sigma_w
real(dp) :: sigma_s
type(Sfilter) :: epsr1,epsr2,ZT
type(Sfilter) :: YC,YD
character(LEN=error_message_length) :: message
! START
! Check the cable parameters
rw=cable%parameters(1)
rd=cable%parameters(2)
s=cable%parameters(3)
rs=cable%parameters(4)
t=cable%parameters(5)
rd2=cable%parameters(6)
sigma_w=cable%parameters(7)
sigma_s=cable%parameters(8)
epsr1=cable%dielectric_filter(1)
epsr2=cable%dielectric_filter(2)
ZT=cable%transfer_impedance(1)
cable_spec_error=.FALSE. ! assume no errors initially
message=''
CALL shielded_twisted_pair_check(rw,rd,s,rs,rd2,cable_spec_error,cable%cable_name,message)
CALL conductivity_check(sigma_w,cable_spec_error,cable%cable_name,message)
CALL conductivity_check(sigma_s,cable_spec_error,cable%cable_name,message)
CALL dielectric_check(epsr1,cable_spec_error,cable%cable_name,message)
CALL dielectric_check(epsr2,cable_spec_error,cable%cable_name,message)
CALL transfer_impedance_check(Zt,cable_spec_error,cable%cable_name,message)
CALL surface_impedance_check(ZT,sigma_s,rs,t,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
! pre-calculate inductance matrix elements for two conductors in a cylindrical shield
domain=1
epsr=evaluate_Sfilter_high_frequency_limit(epsr1)
! First calculate the elements of the 2x2 inductance matrix before
! determining the common mode and differential mode inductance and capacitance
jw=jwA_filter(1d0)
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 shielded twisted pairs'
n_conductors=3
CALL allocate_and_reset_PUL_data(PUL,n_conductors)
PUL%shape(1:n_conductors)=circle
PUL%x(1)=-s/2d0
PUL%y(1)=0.0
PUL%r(1)=rw
PUL%rd(1)=rd
PUL%epsr(1)=epsr1
PUL%x(2)=s/2d0
PUL%y(2)=0.0
PUL%r(2)=rw
PUL%rd(2)=rd
PUL%epsr(2)=epsr1
PUL%epsr_background = 1d0 ! permittivity of homogeneous dielectric medium surrounding conductors (air)
! no ground plane
PUL%ground_plane_present=.FALSE.
! add overshield i.e. the twinax shield
PUL%overshield_present=.TRUE.
PUL%overshield_x = 0d0 ! shield is centred at the origin in this calculation
PUL%overshield_y = 0d0
PUL%overshield_r = rs ! shielded twisted pair shield radius
CALL PUL_LC_Laplace(PUL,cable%cable_name,cable%Y_fit_model_order,cable%Y_fit_freq_spec,domain)
! Theory_Manual_Eqn 3.21
! there may be slight asymmmetry due to meshing so average diagonal and off diagonal elements
L11=(PUL%L%mat(1,1)+PUL%L%mat(1,1))/2d0
L12=(PUL%L%mat(1,2)+PUL%L%mat(2,1))/2d0
C11=(PUL%C%mat(1,1)+PUL%C%mat(1,1))/2d0
C12=(PUL%C%mat(1,2)+PUL%C%mat(2,1))/2d0
dielectric_is_homogeneous=.FALSE.
CALL shielded_twisted_pair_cm_dm_parameter_calculation(L11,L12,C11,C12,epsr,LC,LD,CC,CD,dielectric_is_homogeneous)
! Theory_Manual_Eqn 3.22
YD=0.5d0*( PUL%Yfilter%sfilter_mat(1,1)+((-1d0)*PUL%Yfilter%sfilter_mat(1,2)) )
YC=2.0d0*( PUL%Yfilter%sfilter_mat(1,1)+PUL%Yfilter%sfilter_mat(1,2) )
else
! See C.R. Paul, 1st edition, equation 3.67a,b with cos(thetaij)=-1! Theory_Manual_Eqn 2.27, 2.28
L11=(mu0/(2d0*pi))*log( (rs**2-(s/2d0)**2)/(rs*rw) )
L12=(mu0/(2d0*pi))*log( (s/(2d0*rs)) * (rs**2+(s/2d0)**2)/(2d0*(s/2d0)**2) )
dielectric_is_homogeneous=.TRUE.
CALL shielded_twisted_pair_cm_dm_parameter_calculation(L11,L12,C11,C12,epsr,LC,LD,CC,CD,dielectric_is_homogeneous)
YD=CD*jw
YC=CC*jw
end if
! DOMAIN 1: Set the parameters for the internal differential mode domain
domain=1
cable%n_internal_conductors_in_domain(domain)=2
! The number of modes in the internal differential mode domain is 1
dim=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))
cable%Z_domain(domain)%dim=dim
ALLOCATE(cable%Z_domain(domain)%sfilter_mat(dim,dim))
cable%Y_domain(domain)%dim=dim
ALLOCATE(cable%Y_domain(domain)%sfilter_mat(dim,dim))
cable%L_domain(domain)%mat(1,1)=LD
cable%Z_domain(domain)%sfilter_mat(1,1)=cable%L_domain(domain)%mat(1,1)*jw
cable%C_domain(domain)%mat(1,1)=CD
cable%Y_domain(domain)%sfilter_mat(1,1)=YD
! DOMAIN 2: Set the parameters for the internal common mode domain
domain=2
cable%n_internal_conductors_in_domain(domain)=2
! The number of modes in the internal differential mode domain is 1
dim=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))
cable%Z_domain(domain)%dim=dim
ALLOCATE(cable%Z_domain(domain)%sfilter_mat(dim,dim))
cable%Y_domain(domain)%dim=dim
ALLOCATE(cable%Y_domain(domain)%sfilter_mat(dim,dim))
cable%L_domain(domain)%mat(1,1)=LC
cable%Z_domain(domain)%sfilter_mat(1,1)=cable%L_domain(domain)%mat(1,1)*jw
cable%C_domain(domain)%mat(1,1)=CC
cable%Y_domain(domain)%sfilter_mat(1,1)=YC
if (use_laplace) CALL deallocate_PUL_data(PUL) ! deallocate the PUL data structure
! Set the domain decomposition matrices ! Theory_Manual_Eqn 6.11, 6.12
! The dimension of the domain transformation matrices is 4
dim=4
cable%MI%dim=dim
ALLOCATE(cable%MI%mat(dim,dim))
cable%MV%dim=dim
ALLOCATE(cable%MV%mat(dim,dim))
cable%MI%mat(1,1)=0.5D0
cable%MI%mat(1,2)=-0.5d0
cable%MI%mat(1,3)=0d0
cable%MI%mat(1,4)=0d0
cable%MI%mat(2,1)=1d0
cable%MI%mat(2,2)=1d0
cable%MI%mat(2,3)=0d0
cable%MI%mat(2,4)=0d0
cable%MI%mat(3,1)=1d0
cable%MI%mat(3,2)=1d0
cable%MI%mat(3,3)=1d0
cable%MI%mat(3,4)=0d0
cable%MI%mat(4,1)=1d0
cable%MI%mat(4,2)=1d0
cable%MI%mat(4,3)=1d0
cable%MI%mat(4,4)=1d0
cable%MV%mat(1,1)=1D0
cable%MV%mat(1,2)=-1d0
cable%MV%mat(1,3)=0d0
cable%MV%mat(1,4)=0d0
cable%MV%mat(2,1)=0.5d0
cable%MV%mat(2,2)=0.5d0
cable%MV%mat(2,3)=-1d0
cable%MV%mat(2,4)=0d0
cable%MV%mat(3,1)=0d0
cable%MV%mat(3,2)=0d0
cable%MV%mat(3,3)=1d0
cable%MV%mat(3,4)=-1d0
cable%MV%mat(4,1)=0d0
cable%MV%mat(4,2)=0d0
cable%MV%mat(4,3)=0d0
cable%MV%mat(4,4)=1d0
! Set the local reference conductor numbering
ALLOCATE( cable%local_reference_conductor(3) )
cable%local_reference_conductor(1)=2 ! differential mode, reference is the second conductor
cable%local_reference_conductor(2)=3 ! common mode, reference is the shield conductor
cable%local_reference_conductor(3)=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)=2 ! differential mode domain
cable%local_domain_n_conductors(2)=2 ! common mode
cable%local_domain_n_conductors(3)=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=circle
cable%external_model(1)%conductor_radius=rs
cable%external_model(1)%dielectric_radius=rd2
cable%external_model(1)%dielectric_epsr=epsr2
! set the conductor impedance model for the differential mode
cable%conductor_impedance(1)%impedance_model_type=impedance_model_type_cylindrical_with_conductivity
cable%conductor_impedance(1)%radius=rw
cable%conductor_impedance(1)%conductivity=sigma_w
cable%conductor_impedance(1)%Resistance_multiplication_factor=1.5d0
! set the conductor impedance model for the common mode
cable%conductor_impedance(2)%impedance_model_type=impedance_model_type_cylindrical_with_conductivity
cable%conductor_impedance(2)%radius=rw
cable%conductor_impedance(2)%conductivity=sigma_w
cable%conductor_impedance(2)%Resistance_multiplication_factor=0.5d0
! set the impedance model for the shield conductor
cable%conductor_impedance(3)%impedance_model_type=impedance_model_type_cylindrical_shield
cable%conductor_impedance(3)%radius=rs
cable%conductor_impedance(3)%thickness=t
cable%conductor_impedance(3)%conductivity=sigma_s
cable%conductor_impedance(3)%ZT_filter=ZT
! Deallocate all filters
CALL deallocate_Sfilter(epsr1)
CALL deallocate_Sfilter(epsr2)
CALL deallocate_Sfilter(ZT)
CALL deallocate_Sfilter(jw)
ALLOCATE( cable%conductor_label(1:cable%tot_n_conductors) )
cable%conductor_label(1)='Cable name: '//trim(cable%cable_name)// &
'. type: '//trim(cable%cable_type_string)//'. conductor 1 : Twisted pair wire 1'
cable%conductor_label(2)='Cable name: '//trim(cable%cable_name)// &
'. type: '//trim(cable%cable_type_string)//'. conductor 2 : Twisted pair wire 2'
cable%conductor_label(3)='Cable name: '//trim(cable%cable_name)// &
'. type: '//trim(cable%cable_type_string)//'. conductor 3 : Shield'
END SUBROUTINE shielded_twisted_pair_set_internal_domain_information
!
! NAME
! shielded_twisted_pair_plot
!
! AUTHORS
! Chris Smartt
!
! DESCRIPTION
! plot shielded twisted pair cable
!
! COMMENTS
! the angle has NO impact here due to the twisting
!
! HISTORY
!
! started 5/9/2016 CJS based on shielded_twisted_pair.F90
!
!
SUBROUTINE shielded_twisted_pair_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
real(dp) :: x,y,r,rd
real(dp) :: s
! START
! Inner conductors
r=cable%parameters(1) ! inner conductor radius
rd=cable%parameters(2) ! inner dielecric radius
s=cable%parameters(3) ! inner conductor separation
! plot inner conductor, 1
x=x_offset+(s/2d0)
y=y_offset
CALL write_circle(x,y,r,conductor_geometry_file_unit,xmin,xmax,ymin,ymax)
! plot inner conductor, 2
x=x_offset-(s/2d0)
y=y_offset
CALL write_circle(x,y,r,conductor_geometry_file_unit,xmin,xmax,ymin,ymax)
! plot inner dielectric, 1
x=x_offset+(s/2d0)
y=y_offset
CALL write_circle(x,y,rd,dielectric_geometry_file_unit,xmin,xmax,ymin,ymax)
! plot inner dielectric, 2
x=x_offset-(s/2d0)
y=y_offset
CALL write_circle(x,y,rd,dielectric_geometry_file_unit,xmin,xmax,ymin,ymax)
! plot shield conductor
r=cable%parameters(4) ! outer shield radius
x=x_offset
y=y_offset
CALL write_circle(x,y,r,conductor_geometry_file_unit,xmin,xmax,ymin,ymax)
! plot circular dielectric
rd=cable%parameters(6) ! outer dielectric radius
x=x_offset
y=y_offset
CALL write_circle(x,y,rd,dielectric_geometry_file_unit,xmin,xmax,ymin,ymax)
RETURN
END SUBROUTINE shielded_twisted_pair_plot