Ensamblaje sobre el Vector W
[Ensamblaje los elementos en el sistema lineal]

Diagrama de colaboración para Ensamblaje sobre el Vector W:

Funciones

void bobinaVectorW (double *p_elemento, basicos *p_DATOS)
 Esta Función se encarga del ensamblaje de valores en el vector W de una bobina.
void fuenteVccVectorW (double *p_elemento, basicos *p_DATOS, float t)
 Esta Función se encarga del ensamblaje de valores en el vector W de una fuente independiente de tensión.
void fuenteIVectorW (double *p_elemento, basicos *p_DATOS, float t)
 Esta Función se encarga del ensamblaje de valores en el vector W de una fuente independiente de corriente.

Descripción detallada

Este grupo de funciones es el encargado de realizar el ensamblaje en el Vector W.

Documentación de las funciones

void bobinaVectorW ( double *  p_elemento,
basicos p_DATOS 
)

Esta Función se encarga del ensamblaje de valores en el vector W de una bobina.

Este elemento pertenece al grupo A2 por lo que siempre aporta una ecuación adicional

Parámetros:
elemento Puntero a vector con los valores de definición del elemento
p_DATOS Puntero a estructura basicos

void fuenteIVectorW ( double *  p_elemento,
basicos p_DATOS,
float  t 
)

Esta Función se encarga del ensamblaje de valores en el vector W de una fuente independiente de corriente.

Este elemento pertenece al grupo A2 por lo que siempre aporta una ecuación adicional

Parámetros:
elemento Puntero a vector con los valores de definición del elemento
p_DATOS Puntero a estructura basicos
t Valor de tiempo, sirve para evaluar la excitación en el instante dado

Definición en la línea 2077 del archivo elementos.cpp.

02078 {
02079 float limite_a=p_elemento[10];
02080 float limite_b=p_elemento[11];
02081 double tipo_excitacion0= p_elemento[6];
02082 int tipo_excitacion=int(tipo_excitacion0);
02083 float valor_n,valor_n1,valor;
02084 cout<<"valor de t: "<<t<<endl;
02085 if(p_elemento[2]==1)//Ver si es EA
02086 {
02087 cout<<"Error de concepto en fuente de intensidad "<<p_elemento[0]<<
02088 " con ID relativo: "<<p_elemento[1]<<
02089 ", pertenece al grupo A3 no puede ser ecuación adicional."<<endl;      
02090 }
02091 
02092 switch (tipo_excitacion)
02093 {
02094 
02095         case 0:{if(p_DATOS->transitorio==1){valor_n=p_elemento[7];valor_n1=p_elemento[7];}
02096         if(p_DATOS->estacionario==1){valor=p_elemento[7];}
02097         break;}                 
02098         case 1:{if(p_DATOS->transitorio==1)
02099         {valor_n=p_elemento[12]*sin(p_elemento[9]*2*M_PI*(t*p_DATOS->h)+(p_elemento[8]*RAD));                           
02100         cout<<"valor del seno de valor_n: "<<valor_n<<endl;
02101         valor_n1=p_elemento[12]*sin(p_elemento[9]*2*M_PI*(t*p_DATOS->h+p_DATOS->h));
02102         cout<<"valor del argumento del seno: "<<t*p_DATOS->h+p_DATOS->h<<endl;
02103         cout<<"valor del elemento 8: "<<p_elemento[9]<<endl;
02104         cout<<"valor del seno de valor_n1: "<<valor_n1<<endl;
02105         }
02106         if(p_DATOS->estacionario==1){valor=p_elemento[12];}
02107         break;}//Seno 
02108         case 2:{if(p_DATOS->transitorio==1)
02109         {valor_n=p_elemento[12]*cos(p_elemento[9]*2*M_PI*(t*p_DATOS->h)+(p_elemento[8]*RAD));                           
02110         cout<<"valor del seno: "<<valor_n<<endl;
02111         valor_n1=p_elemento[12]*cos(p_elemento[9]*2*M_PI*(t*p_DATOS->h+p_DATOS->h)+
02112         (p_elemento[8]*RAD));}
02113         if(p_DATOS->estacionario==1){valor=p_elemento[7];}
02114         break;}//Coseno
02115         case 3:{valor_n=p_elemento[12]*t+p_elemento[13];
02116         valor_n1=p_elemento[12]*(t+p_DATOS->h)+p_elemento[13];
02117         break;}//Rampa
02118         case 4:{
02119         if(t<limite_a){valor_n=p_elemento[12];valor_n1=valor_n;}
02120         if(t>=limite_b){valor_n=p_elemento[13];valor_n1=valor_n;}
02121         break;
02122         }//Escalón
02123         case 5:{                
02124         if(t<limite_a){valor_n=p_elemento[14];valor_n1=valor_n;}
02125         if(t>=limite_b){valor_n=p_elemento[15];valor_n1=valor_n;}
02126         if((t>limite_a)&&(t<limite_b))
02127         {
02128         valor_n=p_elemento[12]*t+p_elemento[13];
02129         valor_n1=p_elemento[12]*(t+p_DATOS->h)+p_elemento[13];
02130         }
02131         break;
02132         }//Escalón Modificado
02133 }
02134 
02135 cout<<"valor excitación: "<<valor_n<<endl;
02136         if(p_elemento[4]==0)
02137         {
02138                 if(p_DATOS->transitorio==1){    
02139                 VecSetValue(p_DATOS->Wn,p_elemento[5],valor_n,ADD_VALUES);      
02140                 VecSetValue(p_DATOS->Wn1,p_elemento[5],valor_n1,ADD_VALUES);}
02141                 if(p_DATOS->estacionario==1){   
02142                 VecSetValue(p_DATOS->W,p_elemento[5],valor,ADD_VALUES);}
02143         }
02144         if(p_elemento[5]==0)
02145         {       
02146                 if(p_DATOS->transitorio==1){
02147                 VecSetValue(p_DATOS->Wn,p_elemento[4],-valor_n,ADD_VALUES);
02148                 VecSetValue(p_DATOS->Wn1,p_elemento[4],-valor_n1,ADD_VALUES);}
02149                 if(p_DATOS->estacionario==1){   
02150                 VecSetValue(p_DATOS->W,p_elemento[4],-valor,ADD_VALUES);}       
02151         }
02152         if(p_elemento[4]!=0&&p_elemento[5]!=0)
02153         {
02154                 if(p_DATOS->transitorio==1){
02155                 VecSetValue(p_DATOS->Wn,p_elemento[5],valor_n,ADD_VALUES);      
02156                 VecSetValue(p_DATOS->Wn1,p_elemento[5],valor_n1,ADD_VALUES);
02157                 VecSetValue(p_DATOS->Wn,p_elemento[4],-valor_n,ADD_VALUES);
02158                 VecSetValue(p_DATOS->Wn1,p_elemento[4],-valor_n1,ADD_VALUES);}
02159                 if(p_DATOS->estacionario==1){   
02160                 VecSetValue(p_DATOS->W,p_elemento[5],valor,ADD_VALUES);
02161                 VecSetValue(p_DATOS->W,p_elemento[4],-valor,ADD_VALUES);}
02162         }
02163 }

void fuenteVccVectorW ( double *  p_elemento,
basicos p_DATOS,
float  t 
)

Esta Función se encarga del ensamblaje de valores en el vector W de una fuente independiente de tensión.

Este elemento pertenece al grupo A2 por lo que siempre aporta una ecuación adicional

Parámetros:
elemento Puntero a vector con los valores de definición del elemento
p_DATOS Puntero a estructura basicos
t Valor de tiempo, sirve para evaluar la excitación en el instante dado

Definición en la línea 1995 del archivo elementos.cpp.

01996 {
01997 float limite_a=p_elemento[10];
01998 float limite_b=p_elemento[11];
01999 double tipo_excitacion0= p_elemento[6];
02000 int tipo_excitacion=int(tipo_excitacion0);
02001 float valor_n,valor_n1,valor;
02002 
02003 if(p_elemento[2]!=0)//Ver si no es EA
02004 {
02005 cout<<"Error de concepto en fuente de tensión "<<p_elemento[0]<<
02006 " con ID relativo: "<<p_elemento[1]<<
02007 ", este siempre es EA."<<endl;  
02008 }
02009 
02010 switch (tipo_excitacion)
02011 {
02012 
02013         case 0:{if(p_DATOS->transitorio==1){valor_n=p_elemento[7];valor_n1=p_elemento[7];}
02014         if(p_DATOS->estacionario==1){valor=p_elemento[7];}
02015         break;}                 
02016         case 1:{if(p_DATOS->transitorio==1)
02017         {valor_n=p_elemento[12]*sin(p_elemento[9]*2*M_PI*(t*p_DATOS->h)+(p_elemento[8]*RAD));                           
02018         valor_n1=p_elemento[12]*sin(p_elemento[9]*2*M_PI*(t*p_DATOS->h+p_DATOS->h));                            
02019         }
02020         if(p_DATOS->estacionario==1){valor=p_elemento[12];}
02021         break;}//Seno 
02022         case 2:{if(p_DATOS->transitorio==1)
02023         {valor_n=p_elemento[12]*cos(p_elemento[9]*2*M_PI*(t*p_DATOS->h)+(p_elemento[8]*RAD));                           
02024         valor_n1=p_elemento[12]*cos(p_elemento[9]*2*M_PI*(t*p_DATOS->h+p_DATOS->h)+
02025         (p_elemento[8]*RAD));}
02026         if(p_DATOS->estacionario==1){valor=p_elemento[7];}
02027         break;}//Coseno
02028         case 3:{
02029         valor_n=p_elemento[12]*t+p_elemento[13];
02030         valor_n1=p_elemento[12]*(t+p_DATOS->h)+p_elemento[13];
02031         break;}//Rampa
02032         case 4:{
02033         if(t<limite_a){valor_n=p_elemento[12];valor_n1=valor_n;}
02034         if(t>=limite_a){valor_n=p_elemento[13];valor_n1=valor_n;}
02035         break;
02036         }//Escalón
02037         case 5:{                
02038         if(t<=limite_a){valor_n=p_elemento[14];valor_n1=valor_n;}
02039         if(t>=limite_b){valor_n=p_elemento[15];valor_n1=valor_n;}
02040         if((t>limite_a)&&(t<limite_b))
02041         {
02042         valor_n=p_elemento[12]*t+p_elemento[13];
02043         valor_n1=p_elemento[12]*(t+p_DATOS->h)+p_elemento[13];
02044         }
02045         break;
02046         }//Escalón Modificado
02047         }
02048 
02049         if(p_elemento[4]==0)
02050         {
02051                 if(p_DATOS->transitorio==1){    
02052                 VecSetValue(p_DATOS->Wn,p_DATOS->apuntaA2,valor_n,ADD_VALUES);  
02053                 VecSetValue(p_DATOS->Wn1,p_DATOS->apuntaA2,valor_n1,ADD_VALUES);}
02054                 if(p_DATOS->estacionario==1){   
02055                 VecSetValue(p_DATOS->W,p_DATOS->apuntaA2,valor,ADD_VALUES);}
02056         }
02057         if(p_elemento[5]==0)
02058         {       
02059                 if(p_DATOS->transitorio==1){
02060                 VecSetValue(p_DATOS->Wn,p_DATOS->apuntaA2,valor_n,ADD_VALUES);
02061                 VecSetValue(p_DATOS->Wn1,p_DATOS->apuntaA2,valor_n1,ADD_VALUES);}
02062                 if(p_DATOS->estacionario==1){   
02063                 VecSetValue(p_DATOS->W,p_DATOS->apuntaA2,valor,ADD_VALUES);}    
02064         }
02065         if(p_elemento[4]!=0&&p_elemento[5]!=0)
02066         {
02067                 if(p_DATOS->transitorio==1){
02068                 VecSetValue(p_DATOS->Wn,p_DATOS->apuntaA2,valor_n,ADD_VALUES);
02069                 VecSetValue(p_DATOS->Wn1,p_DATOS->apuntaA2,valor_n1,ADD_VALUES);}
02070                 if(p_DATOS->estacionario==1){   
02071                 VecSetValue(p_DATOS->W,p_DATOS->apuntaA2,valor,ADD_VALUES);}
02072         }
02073 //Incremento en uno el valor para poder poner otro elemento EA si lo hubiese
02074 p_DATOS->apuntaA2=p_DATOS->apuntaA2+1;
02075 }


Generado el Wed Sep 2 16:04:56 2009 para Software de Análisis de Circuitos Lineales mediante Métodos Avanzados de Resolución para Uso Docente. Manual de Código Fuente por  doxygen 1.5.6