/* __h2.c Hopf(neutral saddle) point-->Hopf(neutral saddle) Curve using double-bordering (A^2 + kI) method October 16, 1997 */ #include "common.h" #include "corefunc.h" #include "stagen.h" #include "_util.h" #include "_linalg.h" #include "__h2.h" #include "_conti.h" /* continuer's (i.e. generator's) data */ /* Numbers of visible name classes; the order must be the same as in /names section */ #define CL_TIME 1 #define CL_PHASE 2 #define CL_PAR 3 #define CL_EIGV 4 #define CL_TEST 5 #define CL_UTEST 6 #define _eps staData->eps #define _dx staData->dx #define _branch staData->branch #define PI 4.0*atan(1.0) #define PID2 2.0*atan(1.0) Local(StagenDataPtr) stagenData=NULL; /* global copy of Starter's parameter passed by content */ Local(__h_DataPtr) staData=NULL; /* global copy of our own data area */ Local(StdLinAlgDataPtr) staFunc=NULL; /* global copy of pointer standard linear algebra struct */ Local(UtilitiesDataPtr) staUtil=NULL; /* global copy of pointer to utility functions struct */ Local(Int2Ptr) active=NULL; /* index of active parameter(s) */ Local(Int2) nappend; /* number of appended user functions */ Local(Int2) npar,nphase; /* dims of par and phase spaces */ Local(Vector) x0=NULL; /* approximate first point */ Local(Vector) x1=NULL; /* current point */ Local(Vector) v0=NULL; /* approximate tangent vector at first point */ Local(FloatHiPtr) X=NULL; /* current set of phases used by Decoder */ Local(FloatHiPtr) P=NULL; /* current set of params used by Decoder */ Local(FloatHiPtr) E=NULL; /* current set of eigenvalues */ Local(VAR) zero=0; #include "_symdiff.c" /* Working space for defining and test functions */ Local(Matrix) A,B,C,C1,C2,AA,AT,AH,BB,BP,BT,DE,G,BP1,BT1,JM1,JM2,TVM,TVMT,HM,GM1,GM2,Jac1,Jac2; Local(Vector) vn,xp,xp1,V,W,F,Vp,L,VV,VB,V2,W2,Q1,Q2,Qr,Qi,PP1,PP2,Pr,Pi,aa,bb,cc,dd,rr; Local(Vector) G0,G1,G2,G3,TV1,TV2,HV,HV1,HV2,HW1,HW2,VB1,VB2,V1,W1,RL,K1,K2,K3,OK1,OK2,PIVOT; Local(Vector) Re=NULL,Im=NULL; Local(VAR) test[4],dx,ddx,dddx; Local(Int2) m1,n1,m2,n2; #define BifVector stagenData->BifDataOutPtr #define ActiveUTestFuns staData->ActiveUserTest #define nUTest stagenData->udFuncNum /* Initialize function common for Starter and Decoder */ Local(Int2) Init(StagenDataPtr stagenptr) { Int2 i,m; /* Set local pointers to data area and functions array */ stagenData=stagenptr; staData=(__h_DataPtr)stagenptr->StaPtr; /* Use standard linear algebra */ staFunc=(StdLinAlgDataPtr)stagenptr->StdLinAlg; staUtil=(UtilitiesDataPtr)stagenptr->Utilities; /* Get dimensions of Phase and Par spaces */ npar=*stagenptr->ClassDim[CL_PAR]; nphase=*stagenptr->ClassDim[CL_PHASE]; /* Call udf SetInitPoint function, if any */ if (!stagenData->ContDataGenLen && FuncParActive(staData->SetInitPoint)) { /* there is UDF to set InitPoint */ Int2 Index_X,Index_P; Index_X=stagenptr->ParIndex(staData,X); Index_P=stagenptr->ParIndex(staData,P); FuncParCall(staData->SetInitPoint)(staData); stagenptr->UpdatePar(Index_X); stagenptr->UpdatePar(Index_P); } /* Do all the checks related to active parameters and find them */ for (i=0,nappend=0; iCheckParameters(staData->Active,npar,active,2+nappend)) { active=MemFree(active); return 1; } /* Checks related to analytical/numerical differentiation of RHS */ if (!stagenData->Der1 && staData->dx <= 0) { staUtil->ErrMessage("Analytical Jacoby matrix is undefined.\nRecompile system or set positive 'Increment'"); active=MemFree(active); return 2; } if ((staData->ActiveSing[0] || staData->ActiveSing[1] || staData->ActiveSing[2]) && !stagenData->Der2 && staData->dx <= 0) { staUtil->ErrMessage("Analytical Hessian matrix is undefined.\nRecompile system or set positive 'Increment'"); active=MemFree(active); return 3; } /* Create a copy parameters */ m=nphase*(nphase-1)/2; if (!P) { P=MemNew(sizeof(FloatHi)*npar); MemCopy(P,staData->P,sizeof(FloatHi)*npar); /* Allocate memory for copy of phase vars used by Decoder */ X=MemNew(sizeof(FloatHi)*nphase); /* Allocate memory for eigenvalues */ E=MemNew(sizeof(FloatHi)*2*nphase); Re=staFunc->CreateVector(nphase); Im=staFunc->CreateVector(nphase); /* Allocate memory for working space used in defining and test functions */ A=staFunc->CreateMatrix(nphase,nphase); B=staFunc->CreateMatrix(nphase,nphase); BP=staFunc->CreateMatrix(nphase+2,nphase+2); BT=staFunc->CreateMatrix(nphase+2,nphase+2); BP1=staFunc->CreateMatrix(nphase+1,nphase+1); BT1=staFunc->CreateMatrix(nphase+1,nphase+1); V1=staFunc->CreateVector(nphase+2); W1=staFunc->CreateVector(nphase+2); V2=staFunc->CreateVector(nphase+2); W2=staFunc->CreateVector(nphase+2); VB=staFunc->CreateVector(nphase); VB1=staFunc->CreateVector(nphase+2); VB2=staFunc->CreateVector(nphase+2); VV=staFunc->CreateVector(2*nphase); Q1=staFunc->CreateVector(nphase+2); Q2=staFunc->CreateVector(nphase+2); PP1=staFunc->CreateVector(nphase+2); PP2=staFunc->CreateVector(nphase+2); RL=staFunc->CreateVector(m+1); K1=staFunc->CreateVector(nphase); K2=staFunc->CreateVector(nphase); K3=staFunc->CreateVector(nphase); OK1=staFunc->CreateVector(m+2); OK2=staFunc->CreateVector(m+2); G0=staFunc->CreateVector(nphase+3+nappend); G1=staFunc->CreateVector(nphase+3+nappend); G2=staFunc->CreateVector(nphase+3+nappend); G3=staFunc->CreateVector(nphase+3+nappend); JM1=staFunc->CreateMatrix(nphase+2,nphase+2); JM2=staFunc->CreateMatrix(nphase+3,nphase+3); TV1=staFunc->CreateVector(nphase+2); TV2=staFunc->CreateVector(nphase+3); TVM=staFunc->CreateMatrix(nphase+3,3); TVMT=staFunc->CreateMatrix(3,nphase+3); HM=staFunc->CreateMatrix(3,3); HV=staFunc->CreateVector(3); HV1=staFunc->CreateVector(nphase); HV2=staFunc->CreateVector(nphase); HW1=staFunc->CreateVector(nphase); HW2=staFunc->CreateVector(nphase); GM1=staFunc->CreateMatrix(nphase+3+nappend,4); GM2=staFunc->CreateMatrix(nphase+3+nappend,4); G=staFunc->CreateMatrix(2,2); AA=staFunc->CreateMatrix(nphase,nphase); BB=staFunc->CreateMatrix(2*nphase,2*nphase); C=staFunc->CreateMatrix(nphase+nappend,nphase+3+nappend); C1=staFunc->CreateMatrix(nphase,nphase+3+nappend); C2=staFunc->CreateMatrix(2,nphase+3+nappend); x0=staFunc->CreateVector(nphase+3+nappend); x1=staFunc->CreateVector(nphase+3+nappend); v0=staFunc->CreateVector(nphase+3+nappend); vn=staFunc->CreateVector(nphase); xp=staFunc->CreateVector(nphase+2+nappend); xp1=staFunc->CreateVector(nphase+2+nappend); V=staFunc->CreateVector(nphase); W=staFunc->CreateVector(nphase); F=staFunc->CreateVector(nphase); L=staFunc->CreateVector(nphase); Vp=staFunc->CreateVector(nphase+2+nappend); BifVector=MemNew(sizeof(FloatHi)*(4*nphase+2)); if (nphase==1) staData->ActiveSing[0]=FALSE; if ((nphase>1) && staData->ActiveSing[0]) { /* Bautin detection ON*/ AT=staFunc->CreateMatrix(nphase,nphase); Qr=staFunc->CreateVector(nphase); Qi=staFunc->CreateVector(nphase); Pr=staFunc->CreateVector(nphase); Pi=staFunc->CreateVector(nphase); aa=staFunc->CreateVector(nphase); bb=staFunc->CreateVector(nphase); cc=staFunc->CreateVector(nphase); dd=staFunc->CreateVector(nphase); rr=staFunc->CreateVector(nphase); } if (nphase<=3) staData->ActiveSing[3]=FALSE; if ((nphase>3) && staData->ActiveSing[3]) { /*Double Hopf detection ON*/ AH=staFunc->CreateMatrix(nphase-2,nphase-2); DE=staFunc->CreateMatrix((nphase-2)*(nphase-3)/2,(nphase-2)*(nphase-3)/2); } } Der2num=ind2_(nphase+npar-1,nphase+npar-1)+1; Der3num=ind3_(nphase+npar-1,nphase+npar-1,nphase+npar-1)+1; dx=_dx; ddx=pow(dx,3.0/4.0); dddx=pow(dx,3.0/6.0); ((ContDataPtr)stagenptr->GenPtr)->ActiveSing=staData->ActiveSing; ((ContDataPtr)stagenptr->GenPtr)->ActiveUserTest=staData->ActiveUserTest; ((ContDataPtr)stagenptr->GenPtr)->TypeTestFuns=staData->TypeTestFunc; ((ContDataPtr)stagenptr->GenPtr)->ZeroTestFuns=staData->ZeroTestFunc; return 0; } Local(Int2) JacRHS (Vector x, Matrix J); Local(void) Term (Boolean OK) { Int2 m; stagenData->BifDataOutLen=0; m=nphase*(nphase-1)/2; if (OK&&P&&(stagenData->Reason==RF_CLEAN)) { /* saving global vectors V1,V2,W1 and W2 used in defining functions */ stagenData->ContDataStaLen=4*(nphase+2)*sizeof(FloatHi); stagenData->ContDataStaPtr=MemNew(stagenData->ContDataStaLen); MemCopy(stagenData->ContDataStaPtr,V1,(nphase+2)*sizeof(FloatHi)); MemCopy(stagenData->ContDataStaPtr+(nphase+2),W1,(nphase+2)*sizeof(FloatHi)); MemCopy(stagenData->ContDataStaPtr+2*(nphase+2),V2,(nphase+2)*sizeof(FloatHi)); MemCopy(stagenData->ContDataStaPtr+3*(nphase+2),W2,(nphase+2)*sizeof(FloatHi)); stagenData->BifDataOutLen=sizeof(FloatHi)*(2*nphase+1); MemCopy(BifVector,VB1,sizeof(FloatHi)*nphase); MemCopy(BifVector+nphase,VB2,sizeof(FloatHi)*nphase); BifVector[2*nphase]=x1[nphase+2+nappend]; stagenData->ProcessPoint(stagenData->ProcFuncNum, x1, "Last point"); } else { stagenData->ContDataStaLen=0; stagenData->BifDataOutLen=0; } x0=staFunc->DestroyVector(x0); x1=staFunc->DestroyVector(x1); v0=staFunc->DestroyVector(v0); active=MemFree(active); if (P) { P=MemFree(P); X=MemFree(X); E=MemFree(E); Re=staFunc->DestroyVector(Re); Im=staFunc->DestroyVector(Im); A=staFunc->DestroyMatrix(A); B=staFunc->DestroyMatrix(B); BP=staFunc->DestroyMatrix(BP); BT=staFunc->DestroyMatrix(BT); BP1=staFunc->DestroyMatrix(BP1); BT1=staFunc->DestroyMatrix(BT1); V1=staFunc->DestroyVector(V1); W1=staFunc->DestroyVector(W1); V2=staFunc->DestroyVector(V2); W2=staFunc->DestroyVector(W2); Q1=staFunc->DestroyVector(Q1); Q2=staFunc->DestroyVector(Q2); PP1=staFunc->DestroyVector(PP1); PP2=staFunc->DestroyVector(PP2); AA=staFunc->DestroyMatrix(AA); BB=staFunc->DestroyMatrix(BB); VV=staFunc->DestroyVector(VV); VB=staFunc->DestroyVector(VB); VB1=staFunc->DestroyVector(VB1); VB2=staFunc->DestroyVector(VB2); RL=staFunc->DestroyVector(RL); K1=staFunc->DestroyVector(K1); K2=staFunc->DestroyVector(K2); K3=staFunc->DestroyVector(K3); OK1=staFunc->DestroyVector(OK1); OK2=staFunc->DestroyVector(OK2); G0=staFunc->DestroyVector(G0); G1=staFunc->DestroyVector(G1); G2=staFunc->DestroyVector(G2); G3=staFunc->DestroyVector(G3); JM1=staFunc->DestroyMatrix(JM1); JM2=staFunc->DestroyMatrix(JM2); TV1=staFunc->DestroyVector(TV1); TV2=staFunc->DestroyVector(TV2); TVM=staFunc->DestroyMatrix(TVM); TVMT=staFunc->DestroyMatrix(TVMT); HM=staFunc->DestroyMatrix(HM); HV=staFunc->DestroyVector(HV); HV1=staFunc->DestroyVector(HV1); HV2=staFunc->DestroyVector(HV2); HW1=staFunc->DestroyVector(HW1); HW2=staFunc->DestroyVector(HW2); GM1=staFunc->DestroyMatrix(GM1); GM2=staFunc->DestroyMatrix(GM2); G=staFunc->DestroyMatrix(G); C=staFunc->DestroyMatrix(C); C1=staFunc->DestroyMatrix(C1); C2=staFunc->DestroyMatrix(C2); vn=staFunc->DestroyVector(vn); xp=staFunc->DestroyVector(xp); xp1=staFunc->DestroyVector(xp1); V=staFunc->DestroyVector(V); W=staFunc->DestroyVector(W); F=staFunc->DestroyVector(F); Vp=staFunc->DestroyVector(Vp); L=staFunc->DestroyVector(L); BifVector=MemFree(BifVector); if ((nphase>1) && staData->ActiveSing[0]) { /* Bautin detection ON*/ AT=staFunc->DestroyMatrix(AT); Qr=staFunc->DestroyVector(Qr); Qi=staFunc->DestroyVector(Qi); Pr=staFunc->DestroyVector(Pr); Pi=staFunc->DestroyVector(Pi); aa=staFunc->DestroyVector(aa); bb=staFunc->DestroyVector(bb); cc=staFunc->DestroyVector(cc); dd=staFunc->DestroyVector(dd); rr=staFunc->DestroyVector(rr); } if ((nphase>3) && staData->ActiveSing[3]) { /* Double Hopf detection ON*/ AH=staFunc->DestroyMatrix(AH); DE=staFunc->DestroyMatrix(DE); } } } Local(Int2) JacDefEquilib (Vector x, Matrix J); Entry(Int2) h2_JacDefHopf (Vector x, Matrix J); Entry(Int2) h2_DefHopf (Vector x, Vector J); /************/ /* Starters */ /* h -> h2 Starter */ Entry(Int2) h_h2_Starter(StagenDataPtr stagenptr) { Int2 i,j,k,kp1,l,Ret,p,q; Vector v1,v2; Matrix D,DE,S; VAR t,scal; if (!stagenptr) { /* this is request to terminate */ Term(TRUE); return 0; } if (Init(stagenptr)) return 1; /* initialize local variables */ /* Create the first point for generator */ if (stagenptr->ContDataGenLen) { /* "Extension" */ MemCopy(V1,stagenptr->ContDataStaPtr,sizeof(FloatHi)*(nphase+2)); MemCopy(W1,stagenptr->ContDataStaPtr+(nphase+2),sizeof(FloatHi)*(nphase+2)); MemCopy(V2,stagenptr->ContDataStaPtr+2*(nphase+2),sizeof(FloatHi)*(nphase+2)); MemCopy(W2,stagenptr->ContDataStaPtr+3*(nphase+2),sizeof(FloatHi)*(nphase+2)); MemFree(stagenptr->ContDataStaPtr); stagenptr->ContDataStaLen=0; MemCopy(x0,staData->X,sizeof(FloatHi)*nphase); for (i=0; i<2+nappend; i++) x0[nphase+i]=staData->P[active[i]]; x0[nphase+2+nappend]=BifVector[2*nphase]; } else { v1=staFunc->CreateVector(nphase+3+nappend); v2=staFunc->CreateVector(nphase+3+nappend); D=staFunc->CreateMatrix(nphase+2+nappend,nphase+3+nappend); DE=staFunc->CreateMatrix(nphase+3+nappend,nphase+3+nappend); MemCopy(x0,staData->X,sizeof(FloatHi)*nphase); for (i=0; i<2+nappend; i++) x0[nphase+i]=staData->P[active[i]]; JacRHS(x0,C1); if (stagenData->BifDataInOk) { /* Message(MSG_OK,"Eigenvectors exist"); */ /* bifurcation data exist */ /* BifVector[0] to [nphase-1]: Hopf continuation vector */ /* BifVector[nphase] to [2*nphase-1]: global vector L */ /* BifVector[2*nphase]: -lambda*lambda */ MemCopy(V,stagenData->BifDataInPtr,sizeof(FloatHi)*nphase); MemCopy(L,stagenData->BifDataInPtr+nphase,sizeof(FloatHi)*nphase); x0[nphase+2+nappend]=stagenData->BifDataInPtr[2*nphase]; } else { /* bifurcation data to be computed */ /* Message(MSG_OK,"Eigenvectors have to be computed"); */ Ret=0; {int i,j,imin,jmin; double r,s,d,phi,beta,gamma,lambda,omega; staFunc->CopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); staFunc->EigenValues(AA,Re,Im); /* find the pair of eigenvalues with zero sum */ s=fabs(Re[0]+Re[1])+fabs(Im[0]+Im[1]); imin=0; jmin=1; for (i=0; iCopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); for (i=0; iSngSolve(AA,_eps,V) != nphase-1) { staUtil->ErrMessage("Unable to find (+)eigenvector at the initial point"); Ret=1; goto polufinal; } staFunc->CopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); for (i=0; iSngSolve(AA,_eps,W) != nphase-1) { staUtil->ErrMessage("Unable to find (-)eigenvector at the initial point"); Ret=1; goto polufinal; } staFunc->NormalizeVector(V); staFunc->NormalizeVector(W); staFunc->AddScaledVector(V,W,-1.0,L); staFunc->NormalizeVector(L); staFunc->AddScaledVector(V,W,1.0,V); staFunc->NormalizeVector(V); x0[nphase+2+nappend]=-lambda*lambda; } else { /* Hopf */ lambda=Re[imin]; omega=fabs(Im[imin]); /* compute real and imaginary part of the eigenvector */ staFunc->CopyMatrixBlock(C1,0,0,BB,0,0,nphase,nphase); staFunc->CopyMatrixBlock(C1,0,0,BB,nphase,nphase,nphase,nphase); for (i=0; iSngSolve(BB,_eps,VV) != 2*nphase-2) { staUtil->ErrMessage("Unable to find eigenvectors at the initial point"); Ret=1; goto polufinal; } staFunc->CopyVectorElements(VV,0,V,0,nphase); staFunc->CopyVectorElements(VV,nphase,W,0,nphase); /* normalize real and imaginary parts */ d=staFunc->ScalProd(V,V); s=staFunc->ScalProd(W,W); r=staFunc->ScalProd(V,W); if (r!=0) { beta=sqrt(fabs(4.0*r*r+(s-d)*(s-d))); if (2.0*fabs(2.0*r) > beta) { phi=PID2-atan(fabs((s-d)/(2.0*r))); } else { phi=atan(fabs(2.0*r/(s-d))); } if (r < 0) phi=PI-phi; if ((s-d) < 0) phi=-phi; phi/=2.0; for (i=0; iCopyVectorElements(VV,0,V,0,nphase); d=staFunc->ScalProd(V,V); gamma=1/sqrt(d); staFunc->MultiplyVectorScalar(VV,gamma,VV); staFunc->CopyVectorElements(VV,0,V,0,nphase); staFunc->CopyVectorElements(VV,nphase,L,0,nphase); x0[nphase+2+nappend]=omega*omega; staFunc->NormalizeVector(L); } polufinal: if (Ret) { v1=staFunc->DestroyVector(v1); v2=staFunc->DestroyVector(v2); D=staFunc->DestroyMatrix(D); DE=staFunc->DestroyMatrix(DE); Term(FALSE); return 1; } } } if (JacRHS(x0,C1)) return 1; /* for safety */ staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->MultiplyMatrixMatrix(A,A,BP); for (i=0;i= 2) { staFunc->ClearMatrix(BP1); staFunc->CopyMatrixBlock(BP,0,0,BP1,0,0,nphase,nphase); for (i=0; i t) { p = i; q = j; t = fabs(BP1[i][j]); } } } if (k!=q) { /* switch colom k with colom q */ for (i=0;i<=nphase;i++) { t=BP1[i][k]; BP1[i][k]=BP1[i][q]; BP1[i][q]=t; } } if (k!=p) { /* switch row k with row p */ for (i=0;i<=nphase;i++) { t=BP1[k][i]; BP1[k][i]=BP1[p][i]; BP1[p][i]=t; } } if (fabs(BP1[k][k]) == 0.0) { staUtil->ErrMessage("(A^2 + k I) has rank defect at least three"); Term(FALSE); return 1; } kp1=k+1; for (i=kp1;iClearVector(V1); staFunc->ClearVector(V2); staFunc->ClearVector(W1); staFunc->ClearVector(W2); staFunc->CopyMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(VB1); VB1[nphase]=1.0; staFunc->Solve(BT,VB1); staFunc->CopyVector(VB1,V1); staFunc->ClearVector(VB2); VB2[nphase+1]=1.0; staFunc->Solve(BT,VB2); staFunc->CopyVector(VB2,V2); /* solve transposed systems to obtain new values for W1 and W2*/ staFunc->TransposeMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(PP1); PP1[nphase]=1.0; staFunc->Solve(BT,PP1); staFunc->CopyVector(PP1,W1); staFunc->ClearVector(PP2); PP2[nphase+1]=1.0; staFunc->Solve(BT,PP2); staFunc->CopyVector(PP2,W2); for (i=0;iClearVector(V1); staFunc->ClearVector(V2); staFunc->ClearVector(W1); staFunc->ClearVector(W2); V1[0]=1.0; W1[0]=1.0; W2[nphase+1]=1.0; } staFunc->CopyMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(VB1); VB1[nphase]=1.0; staFunc->Solve(BT,VB1); staFunc->ClearVector(VB2); VB2[nphase+1]=1.0; staFunc->Solve(BT,VB2); staFunc->TransposeMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(PP1); PP1[nphase]=1.0; staFunc->Solve(BT,PP1); staFunc->ClearVector(PP2); PP2[nphase+1]=1.0; staFunc->Solve(BT,PP2); /* compute derivatives of g11,g12,g21,g22 */ staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); if (stagenData->Der2) { CopyToX(x0); CopyToP(x0); H=stagenData->Der2(X,P,&zero); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->MultiplyMatrixVector(A,HV1,HV1); staFunc->MultiplyMatrixVector(A,HV2,HV2); staFunc->TransposeMatrix(A,B); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); staFunc->MultiplyMatrixVector(B,HW1,HW1); staFunc->MultiplyMatrixVector(B,HW2,HW2); S=staFunc->CreateMatrix(2,2); for (l=0; lDestroyMatrix(S); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); G0[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV1); G1[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV2); G2[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV1); G3[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV2); } else { staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->MultiplyMatrixVector(A,HV1,HV1); staFunc->MultiplyMatrixVector(A,HV2,HV2); staFunc->TransposeMatrix(A,B); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); staFunc->MultiplyMatrixVector(B,HW1,HW1); staFunc->MultiplyMatrixVector(B,HW2,HW2); Jac1=staFunc->CreateMatrix(nphase,nphase); Jac2=staFunc->CreateMatrix(nphase,nphase); for (l=0; lClearVector(xp); staFunc->CopyVectorElements(x0,0,xp,0,nphase+2); xp[l]+=_dx; staFunc->ClearVector(xp1); staFunc->CopyVectorElements(x0,0,xp1,0,nphase+2); xp1[l]-=_dx; JacRHS(xp,C1); staFunc->CopyMatrixBlock(C1,0,0,Jac1,0,0,nphase,nphase); JacRHS(xp1,C1); staFunc->CopyMatrixBlock(C1,0,0,Jac2,0,0,nphase,nphase); /* g11 */ staFunc->CopyVectorElements(VB1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G0[l]=-staFunc->ScalProd(HW1,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G0[l]-=staFunc->ScalProd(K3,HV1); /* g12 */ staFunc->CopyVectorElements(VB2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G1[l]=-staFunc->ScalProd(HW1,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G1[l]-=staFunc->ScalProd(K3,HV2); /* g21 */ staFunc->CopyVectorElements(VB1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G2[l]=-staFunc->ScalProd(HW2,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G2[l]-=staFunc->ScalProd(K3,HV1); /* g22 */ staFunc->CopyVectorElements(VB2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G3[l]=-staFunc->ScalProd(HW2,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G3[l]-=staFunc->ScalProd(K3,HV2); } Jac1=staFunc->DestroyMatrix(Jac1); Jac2=staFunc->DestroyMatrix(Jac2); G0[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV1); G1[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV2); G2[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV1); G3[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV2); } staFunc->ClearMatrix(BT); staFunc->ClearMatrix(JM1); staFunc->CopyMatrixBlock(C1,0,0,BT,0,0,nphase,nphase+2); /* placement of ones determined by partial pivoting */ PIVOT=staFunc->CreateVector(nphase+2); staFunc->ClearVector(PIVOT); staFunc->TransposeMatrix(BT,JM1); PIVOT[nphase+1]=1.0; for (k=0;k<=nphase;k++) { kp1=k+1; p = k; for (i=kp1;i<=nphase+1;i++) { if (fabs(JM1[i][k]) > fabs(JM1[p][k])) p = i; } PIVOT[k]=p; if (p!=k) PIVOT[nphase+1]=-PIVOT[nphase+1]; t=JM1[p][k]; JM1[p][k]=JM1[k][k]; JM1[k][k]=t; if (t==0.0) continue; for (i=kp1;i<=nphase+1;i++) JM1[i][k]=-JM1[i][k]/t; for (j=kp1;j<=nphase+1;j++) { t=JM1[p][j]; JM1[p][j]=JM1[k][j]; JM1[k][j]=t; if (t!=0.0) for (i=kp1;i<=nphase+1;i++) JM1[i][j]+=JM1[i][k]*t; } } staFunc->ClearVector(TV1); for (i=0;iDestroyVector(PIVOT); /* end of pivot information */ staFunc->ClearMatrix(JM1); staFunc->CopyMatrixBlock(C1,0,0,JM1,0,0,nphase,nphase+2); JM1[nphase][p]=1.0; JM1[nphase+1][q]=1.0; staFunc->ClearVector(TV1); TV1[nphase]=1.0; staFunc->Decomp(JM1); staFunc->Solve(JM1,TV1); staFunc->ClearMatrix(TVM); for (i=0;iClearVector(TV1); TV1[nphase+1]=1.0; staFunc->Solve(JM1,TV1); for (i=0;iTransposeMatrix(TVM,TVMT); staFunc->MultiplyMatrixMatrix(TVMT,TVM,HM); staFunc->Decomp(HM); staFunc->MultiplyMatrixVector(TVMT,G0,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G0); staFunc->MultiplyMatrixVector(TVMT,G1,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G1); staFunc->MultiplyMatrixVector(TVMT,G2,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G2); staFunc->MultiplyMatrixVector(TVMT,G3,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G3); /* Determine component with maximum norm => gives m1,n1 */ G[0][0] = staFunc->NormOfVector(G0); G[0][1] = staFunc->NormOfVector(G1); G[1][0] = staFunc->NormOfVector(G2); G[1][1] = staFunc->NormOfVector(G3); m1=0; n1=0; if (G[0][1] > G[m1][n1]) n1=1; if (G[1][0] > G[m1][n1]) { m1=1; n1=0; } if (G[1][1] > G[m1][n1]) { m1=1; n1=1; } for (i=0;iNormalizeVector(G0); /* G1,G2,G3 are projected on G0 */ for (j=0;j<4;j++) { if (j != (2*m1+n1)) { for (i=0;iScalProd(G1,G0); staFunc->AddScaledVector(G1,G0,-scal,G1); for (i=0;i gives m2,n2 */ for (i=0;iNormOfVector(G0); G[0][1] = staFunc->NormOfVector(G1); G[1][0] = staFunc->NormOfVector(G2); G[1][1] = staFunc->NormOfVector(G3); G[m1][n1]=-1.0; m2=0; n2=0; if (G[0][1] > G[m2][n2]) n2=1; if (G[1][0] > G[m2][n2]) { m2=1; n2=0; } if (G[1][1] > G[m2][n2]) { m2=1; n2=1; } /* Compute a tangent vector to the Hopf curve */ h2_DefHopf(x0,v0); Ret=h2_JacDefHopf(x0,D); for (i=0; iAppendMatrixVector(D,v1,DE); Ret=staFunc->Decomp(DE); if (Ret == 0) { v2[nphase+2+nappend]=stagenData->Forward ? 1 : -1; staFunc->Solve(DE,v2); staFunc->CopyVector(v2,v0); break; }; v1[i]=0.0; }; if (Ret) { staUtil->ErrMessage ("Unable to find a tangent vector to the Hopf curve"); v1=staFunc->DestroyVector(v1); v2=staFunc->DestroyVector(v2); D=staFunc->DestroyMatrix(D); DE=staFunc->DestroyMatrix(DE); Term(FALSE); return 1; } staFunc->NormalizeVector(v0); v1=staFunc->DestroyVector(v1); v2=staFunc->DestroyVector(v2); D=staFunc->DestroyMatrix(D); DE=staFunc->DestroyMatrix(DE); } /* Initialize continuation */ ((ContDataPtr)stagenptr->GenPtr)->X0=x0; ((ContDataPtr)stagenptr->GenPtr)->V0=v0; ((ContDataPtr)stagenptr->GenPtr)->ActiveSing=staData->ActiveSing; ((ContDataPtr)stagenptr->GenPtr)->ActiveUserTest=staData->ActiveUserTest; ((ContDataPtr)stagenptr->GenPtr)->TypeTestFuns=staData->TypeTestFunc; ((ContDataPtr)stagenptr->GenPtr)->ZeroTestFuns=staData->ZeroTestFunc; /* That's all */ return 0; } /********************/ /* zh -> h2 Starter */ Entry(Int2) zh_h2_Starter(StagenDataPtr stagenptr) { Int2 i,j,k,kp1,l,Ret,p,q; Vector v1,v2; Matrix D,DE,S; VAR t,scal; if (!stagenptr) { /* this is request to terminate */ Term(TRUE); return 0; } if (Init(stagenptr)) return 1; /* initialize local variables */ /* Create the first point for generator */ if (stagenptr->ContDataGenLen) { /* "Extension" */ MemCopy(V1,stagenptr->ContDataStaPtr,sizeof(FloatHi)*(nphase+2)); MemCopy(W1,stagenptr->ContDataStaPtr+(nphase+2),sizeof(FloatHi)*(nphase+2)); MemCopy(V2,stagenptr->ContDataStaPtr+2*(nphase+2),sizeof(FloatHi)*(nphase+2)); MemCopy(W2,stagenptr->ContDataStaPtr+3*(nphase+2),sizeof(FloatHi)*(nphase+2)); MemFree(stagenptr->ContDataStaPtr); stagenptr->ContDataStaLen=0; MemCopy(x0,staData->X,sizeof(FloatHi)*nphase); for (i=0; i<2+nappend; i++) x0[nphase+i]=staData->P[active[i]]; x0[nphase+2+nappend]=BifVector[2*nphase]; } else { v1=staFunc->CreateVector(nphase+3+nappend); v2=staFunc->CreateVector(nphase+3+nappend); D=staFunc->CreateMatrix(nphase+2+nappend,nphase+3+nappend); DE=staFunc->CreateMatrix(nphase+3+nappend,nphase+3+nappend); MemCopy(x0,staData->X,sizeof(FloatHi)*nphase); for (i=0; i<2+nappend; i++) x0[nphase+i]=staData->P[active[i]]; JacRHS(x0,C1); if (stagenData->BifDataInOk) { /* Message(MSG_OK,"Eigenvectors exist"); */ /* bifurcation data exist */ /* BifVector[0] to [nphase-1]: Hopf continuation vector */ /* BifVector[nphase] to [2*nphase-1]: global vector L */ /* BifVector[2*nphase]: -lambda*lambda */ MemCopy(V,stagenData->BifDataInPtr,sizeof(FloatHi)*nphase); MemCopy(L,stagenData->BifDataInPtr+nphase,sizeof(FloatHi)*nphase); x0[nphase+2+nappend]=stagenData->BifDataInPtr[2*nphase]; } else { /* bifurcation data to be computed */ /* Message(MSG_OK,"Eigenvectors have to be computed"); */ Ret=0; {int i,j,imin,jmin; double r,s,d,phi,beta,gamma,lambda,omega; staFunc->CopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); staFunc->EigenValues(AA,Re,Im); /* find the pair of eigenvalues with zero sum */ s=fabs(Re[0]+Re[1])+fabs(Im[0]+Im[1]); imin=0; jmin=1; for (i=0; iCopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); for (i=0; iSngSolve(AA,_eps,V) != nphase-1) { staUtil->ErrMessage("Unable to find (+)eigenvector at the initial point"); Ret=1; goto polufinal; } staFunc->CopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); for (i=0; iSngSolve(AA,_eps,W) != nphase-1) { staUtil->ErrMessage("Unable to find (-)eigenvector at the initial point"); Ret=1; goto polufinal; } staFunc->NormalizeVector(V); staFunc->NormalizeVector(W); staFunc->AddScaledVector(V,W,-1.0,L); staFunc->NormalizeVector(L); staFunc->AddScaledVector(V,W,1.0,V); staFunc->NormalizeVector(V); x0[nphase+2+nappend]=-lambda*lambda; } else { /* Hopf */ lambda=Re[imin]; omega=fabs(Im[imin]); /* compute real and imaginary part of the eigenvector */ staFunc->CopyMatrixBlock(C1,0,0,BB,0,0,nphase,nphase); staFunc->CopyMatrixBlock(C1,0,0,BB,nphase,nphase,nphase,nphase); for (i=0; iSngSolve(BB,_eps,VV) != 2*nphase-2) { staUtil->ErrMessage("Unable to find eigenvectors at the initial point"); Ret=1; goto polufinal; } staFunc->CopyVectorElements(VV,0,V,0,nphase); staFunc->CopyVectorElements(VV,nphase,W,0,nphase); /* normalize real and imaginary parts */ d=staFunc->ScalProd(V,V); s=staFunc->ScalProd(W,W); r=staFunc->ScalProd(V,W); if (r!=0) { beta=sqrt(fabs(4.0*r*r+(s-d)*(s-d))); if (2.0*fabs(2.0*r) > beta) { phi=PID2-atan(fabs((s-d)/(2.0*r))); } else { phi=atan(fabs(2.0*r/(s-d))); } if (r < 0) phi=PI-phi; if ((s-d) < 0) phi=-phi; phi/=2.0; for (i=0; iCopyVectorElements(VV,0,V,0,nphase); d=staFunc->ScalProd(V,V); gamma=1/sqrt(d); staFunc->MultiplyVectorScalar(VV,gamma,VV); staFunc->CopyVectorElements(VV,0,V,0,nphase); staFunc->CopyVectorElements(VV,nphase,L,0,nphase); x0[nphase+2+nappend]=omega*omega; staFunc->NormalizeVector(L); } polufinal: if (Ret) { v1=staFunc->DestroyVector(v1); v2=staFunc->DestroyVector(v2); D=staFunc->DestroyMatrix(D); DE=staFunc->DestroyMatrix(DE); Term(FALSE); return 1; } } } if (JacRHS(x0,C1)) return 1; /* for safety */ staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->ClearMatrix(BP); staFunc->MultiplyMatrixMatrix(A,A,BP); for (i=0;i= 2) { staFunc->ClearMatrix(BP1); staFunc->CopyMatrixBlock(BP,0,0,BP1,0,0,nphase,nphase); for (i=0; i t) { p = i; q = j; t = fabs(BP1[i][j]); } } } if (k!=q) { /* switch colom k with colom q */ for (i=0;i<=nphase;i++) { t=BP1[i][k]; BP1[i][k]=BP1[i][q]; BP1[i][q]=t; } } if (k!=p) { /* switch row k with row p */ for (i=0;i<=nphase;i++) { t=BP1[k][i]; BP1[k][i]=BP1[p][i]; BP1[p][i]=t; } } if (fabs(BP1[k][k]) == 0.0) { staUtil->ErrMessage("(A^2 + k I) has rank defect at least three"); Term(FALSE); return 1; } kp1=k+1; for (i=kp1;iClearVector(V1); staFunc->ClearVector(V2); staFunc->ClearVector(W1); staFunc->ClearVector(W2); staFunc->CopyMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(VB1); VB1[nphase]=1.0; staFunc->Solve(BT,VB1); staFunc->CopyVector(VB1,V1); staFunc->ClearVector(VB2); VB2[nphase+1]=1.0; staFunc->Solve(BT,VB2); staFunc->CopyVector(VB2,V2); /* solve transposed systems to obtain new values for W1 and W2*/ staFunc->TransposeMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(PP1); PP1[nphase]=1.0; staFunc->Solve(BT,PP1); staFunc->CopyVector(PP1,W1); staFunc->ClearVector(PP2); PP2[nphase+1]=1.0; staFunc->Solve(BT,PP2); staFunc->CopyVector(PP2,W2); for (i=0;iClearVector(V1); staFunc->ClearVector(V2); staFunc->ClearVector(W1); staFunc->ClearVector(W2); V1[0]=1.0; W1[0]=1.0; W2[nphase+1]=1.0; } staFunc->CopyMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(VB1); VB1[nphase]=1.0; staFunc->Solve(BT,VB1); staFunc->ClearVector(VB2); VB2[nphase+1]=1.0; staFunc->Solve(BT,VB2); staFunc->TransposeMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(PP1); PP1[nphase]=1.0; staFunc->Solve(BT,PP1); staFunc->ClearVector(PP2); PP2[nphase+1]=1.0; staFunc->Solve(BT,PP2); /* compute derivatives of g11,g12,g21,g22 */ staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); if (stagenData->Der2) { CopyToX(x0); CopyToP(x0); H=stagenData->Der2(X,P,&zero); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->MultiplyMatrixVector(A,HV1,HV1); staFunc->MultiplyMatrixVector(A,HV2,HV2); staFunc->TransposeMatrix(A,B); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); staFunc->MultiplyMatrixVector(B,HW1,HW1); staFunc->MultiplyMatrixVector(B,HW2,HW2); S=staFunc->CreateMatrix(2,2); for (l=0; lDestroyMatrix(S); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); G0[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV1); G1[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV2); G2[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV1); G3[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV2); } else { staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->MultiplyMatrixVector(A,HV1,HV1); staFunc->MultiplyMatrixVector(A,HV2,HV2); staFunc->TransposeMatrix(A,B); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); staFunc->MultiplyMatrixVector(B,HW1,HW1); staFunc->MultiplyMatrixVector(B,HW2,HW2); Jac1=staFunc->CreateMatrix(nphase,nphase); Jac2=staFunc->CreateMatrix(nphase,nphase); for (l=0; lClearVector(xp); staFunc->CopyVectorElements(x0,0,xp,0,nphase+2); xp[l]+=_dx; staFunc->ClearVector(xp1); staFunc->CopyVectorElements(x0,0,xp1,0,nphase+2); xp1[l]-=_dx; JacRHS(xp,C1); staFunc->CopyMatrixBlock(C1,0,0,Jac1,0,0,nphase,nphase); JacRHS(xp1,C1); staFunc->CopyMatrixBlock(C1,0,0,Jac2,0,0,nphase,nphase); /* g11 */ staFunc->CopyVectorElements(VB1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G0[l]=-staFunc->ScalProd(HW1,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G0[l]-=staFunc->ScalProd(K3,HV1); /* g12 */ staFunc->CopyVectorElements(VB2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G1[l]=-staFunc->ScalProd(HW1,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G1[l]-=staFunc->ScalProd(K3,HV2); /* g21 */ staFunc->CopyVectorElements(VB1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G2[l]=-staFunc->ScalProd(HW2,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G2[l]-=staFunc->ScalProd(K3,HV1); /* g22 */ staFunc->CopyVectorElements(VB2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G3[l]=-staFunc->ScalProd(HW2,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G3[l]-=staFunc->ScalProd(K3,HV2); } Jac1=staFunc->DestroyMatrix(Jac1); Jac2=staFunc->DestroyMatrix(Jac2); G0[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV1); G1[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV2); G2[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV1); G3[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV2); } staFunc->ClearMatrix(BT); staFunc->ClearMatrix(JM1); staFunc->CopyMatrixBlock(C1,0,0,BT,0,0,nphase,nphase+2); /* placement of ones determined by partial pivoting */ PIVOT=staFunc->CreateVector(nphase+2); staFunc->ClearVector(PIVOT); staFunc->TransposeMatrix(BT,JM1); PIVOT[nphase+1]=1.0; for (k=0;k<=nphase;k++) { kp1=k+1; p = k; for (i=kp1;i<=nphase+1;i++) { if (fabs(JM1[i][k]) > fabs(JM1[p][k])) p = i; } PIVOT[k]=p; if (p!=k) PIVOT[nphase+1]=-PIVOT[nphase+1]; t=JM1[p][k]; JM1[p][k]=JM1[k][k]; JM1[k][k]=t; if (t==0.0) continue; for (i=kp1;i<=nphase+1;i++) JM1[i][k]=-JM1[i][k]/t; for (j=kp1;j<=nphase+1;j++) { t=JM1[p][j]; JM1[p][j]=JM1[k][j]; JM1[k][j]=t; if (t!=0.0) for (i=kp1;i<=nphase+1;i++) JM1[i][j]+=JM1[i][k]*t; } } staFunc->ClearVector(TV1); for (i=0;iDestroyVector(PIVOT); /* end of pivot information */ staFunc->ClearMatrix(JM1); staFunc->CopyMatrixBlock(C1,0,0,JM1,0,0,nphase,nphase+2); JM1[nphase][p]=1.0; JM1[nphase+1][q]=1.0; staFunc->ClearVector(TV1); TV1[nphase]=1.0; Ret=staFunc->Decomp(JM1); if (Ret) return 1; staFunc->Solve(JM1,TV1); staFunc->ClearMatrix(TVM); for (i=0;iClearVector(TV1); TV1[nphase+1]=1.0; staFunc->Solve(JM1,TV1); for (i=0;iTransposeMatrix(TVM,TVMT); staFunc->MultiplyMatrixMatrix(TVMT,TVM,HM); Ret=staFunc->Decomp(HM); if (Ret) return 1; staFunc->MultiplyMatrixVector(TVMT,G0,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G0); staFunc->MultiplyMatrixVector(TVMT,G1,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G1); staFunc->MultiplyMatrixVector(TVMT,G2,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G2); staFunc->MultiplyMatrixVector(TVMT,G3,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G3); /* Determine component with maximum norm => gives m1,n1 */ G[0][0] = staFunc->NormOfVector(G0); G[0][1] = staFunc->NormOfVector(G1); G[1][0] = staFunc->NormOfVector(G2); G[1][1] = staFunc->NormOfVector(G3); m1=0; n1=0; if (G[0][1] > G[m1][n1]) n1=1; if (G[1][0] > G[m1][n1]) { m1=1; n1=0; } if (G[1][1] > G[m1][n1]) { m1=1; n1=1; } for (i=0;iNormalizeVector(G0); /* G1,G2,G3 are projected on G0 */ for (j=0;j<4;j++) { if (j != (2*m1+n1)) { for (i=0;iScalProd(G1,G0); staFunc->AddScaledVector(G1,G0,-scal,G1); for (i=0;i gives m2,n2 */ for (i=0;iNormOfVector(G0); G[0][1] = staFunc->NormOfVector(G1); G[1][0] = staFunc->NormOfVector(G2); G[1][1] = staFunc->NormOfVector(G3); G[m1][n1]=-1.0; m2=0; n2=0; if (G[0][1] > G[m2][n2]) n2=1; if (G[1][0] > G[m2][n2]) { m2=1; n2=0; } if (G[1][1] > G[m2][n2]) { m2=1; n2=1; } /* Compute a tangent vector to the Hopf curve */ h2_DefHopf(x0,v0); Ret=h2_JacDefHopf(x0,D); for (i=0; iAppendMatrixVector(D,v1,DE); Ret=staFunc->Decomp(DE); if (Ret == 0) { v2[nphase+2+nappend]=stagenData->Forward ? 1 : -1; staFunc->Solve(DE,v2); staFunc->CopyVector(v2,v0); break; }; v1[i]=0.0; }; if (Ret) { staUtil->ErrMessage ("Unable to find a tangent vector to the Hopf curve"); v1=staFunc->DestroyVector(v1); v2=staFunc->DestroyVector(v2); D=staFunc->DestroyMatrix(D); DE=staFunc->DestroyMatrix(DE); Term(FALSE); return 1; } staFunc->NormalizeVector(v0); v1=staFunc->DestroyVector(v1); v2=staFunc->DestroyVector(v2); D=staFunc->DestroyMatrix(D); DE=staFunc->DestroyMatrix(DE); } /* Initialize continuation */ ((ContDataPtr)stagenptr->GenPtr)->X0=x0; ((ContDataPtr)stagenptr->GenPtr)->V0=v0; ((ContDataPtr)stagenptr->GenPtr)->ActiveSing=staData->ActiveSing; ((ContDataPtr)stagenptr->GenPtr)->ActiveUserTest=staData->ActiveUserTest; ((ContDataPtr)stagenptr->GenPtr)->TypeTestFuns=staData->TypeTestFunc; ((ContDataPtr)stagenptr->GenPtr)->ZeroTestFuns=staData->ZeroTestFunc; /* That's all */ return 0; } /********************/ /* bt -> h2 Starter */ Entry(Int2) bt_h2_Starter(StagenDataPtr stagenptr) { Int2 i,j,k,kp1,l,Ret,p,q; Vector v1,v2; Matrix D,DE,S; VAR t,scal; if (!stagenptr) { /* this is request to terminate */ Term(TRUE); return 0; } if (Init(stagenptr)) return 1; /* initialize local variables */ /* Create the first point for generator */ if (stagenptr->ContDataGenLen) { /* "Extension" */ MemCopy(V1,stagenptr->ContDataStaPtr,sizeof(FloatHi)*(nphase+2)); MemCopy(W1,stagenptr->ContDataStaPtr+(nphase+2),sizeof(FloatHi)*(nphase+2)); MemCopy(V2,stagenptr->ContDataStaPtr+2*(nphase+2),sizeof(FloatHi)*(nphase+2)); MemCopy(W2,stagenptr->ContDataStaPtr+3*(nphase+2),sizeof(FloatHi)*(nphase+2)); MemFree(stagenptr->ContDataStaPtr); stagenptr->ContDataStaLen=0; MemCopy(x0,staData->X,sizeof(FloatHi)*nphase); for (i=0; i<2+nappend; i++) x0[nphase+i]=staData->P[active[i]]; x0[nphase+2+nappend]=BifVector[2*nphase]; } else { v1=staFunc->CreateVector(nphase+3+nappend); v2=staFunc->CreateVector(nphase+3+nappend); D=staFunc->CreateMatrix(nphase+2+nappend,nphase+3+nappend); DE=staFunc->CreateMatrix(nphase+3+nappend,nphase+3+nappend); MemCopy(x0,staData->X,sizeof(FloatHi)*nphase); for (i=0; i<2+nappend; i++) x0[nphase+i]=staData->P[active[i]]; JacRHS(x0,C1); if (stagenData->BifDataInOk) { /* Message(MSG_OK,"Eigenvectors exist"); */ /* bifurcation data exist */ /* BifVector[0] to [nphase-1]: Hopf continuation vector */ /* BifVector[nphase] to [2*nphase-1]: global vector L */ /* BifVector[2*nphase]: -lambda*lambda */ MemCopy(V,stagenData->BifDataInPtr,sizeof(FloatHi)*nphase); MemCopy(L,stagenData->BifDataInPtr+nphase,sizeof(FloatHi)*nphase); x0[nphase+2+nappend]=stagenData->BifDataInPtr[2*nphase]; } else { /* bifurcation data to be computed */ /* Message(MSG_OK,"Eigenvectors have to be computed"); */ Ret=0; {int i,j,imin,jmin; double r,s,d,phi,beta,gamma,lambda,omega; staFunc->CopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); staFunc->EigenValues(AA,Re,Im); /* find the pair of eigenvalues with zero sum */ s=fabs(Re[0]+Re[1])+fabs(Im[0]+Im[1]); imin=0; jmin=1; for (i=0; iCopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); for (i=0; iSngSolve(AA,_eps,V) != nphase-1) { staUtil->ErrMessage("Unable to find (+)eigenvector at the initial point"); Ret=1; goto polufinal; } staFunc->CopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); for (i=0; iSngSolve(AA,_eps,W) != nphase-1) { staUtil->ErrMessage("Unable to find (-)eigenvector at the initial point"); Ret=1; goto polufinal; } staFunc->NormalizeVector(V); staFunc->NormalizeVector(W); staFunc->AddScaledVector(V,W,-1.0,L); staFunc->NormalizeVector(L); staFunc->AddScaledVector(V,W,1.0,V); staFunc->NormalizeVector(V); x0[nphase+2+nappend]=-lambda*lambda; } else { /* Hopf */ lambda=Re[imin]; omega=fabs(Im[imin]); /* compute real and imaginary part of the eigenvector */ staFunc->CopyMatrixBlock(C1,0,0,BB,0,0,nphase,nphase); staFunc->CopyMatrixBlock(C1,0,0,BB,nphase,nphase,nphase,nphase); for (i=0; iSngSolve(BB,_eps,VV) != 2*nphase-2) { staUtil->ErrMessage("Unable to find eigenvectors at the initial point"); Ret=1; goto polufinal; } staFunc->CopyVectorElements(VV,0,V,0,nphase); staFunc->CopyVectorElements(VV,nphase,W,0,nphase); /* normalize real and imaginary parts */ d=staFunc->ScalProd(V,V); s=staFunc->ScalProd(W,W); r=staFunc->ScalProd(V,W); if (r!=0) { beta=sqrt(fabs(4.0*r*r+(s-d)*(s-d))); if (2.0*fabs(2.0*r) > beta) { phi=PID2-atan(fabs((s-d)/(2.0*r))); } else { phi=atan(fabs(2.0*r/(s-d))); } if (r < 0) phi=PI-phi; if ((s-d) < 0) phi=-phi; phi/=2.0; for (i=0; iCopyVectorElements(VV,0,V,0,nphase); d=staFunc->ScalProd(V,V); gamma=1/sqrt(d); staFunc->MultiplyVectorScalar(VV,gamma,VV); staFunc->CopyVectorElements(VV,0,V,0,nphase); staFunc->CopyVectorElements(VV,nphase,L,0,nphase); x0[nphase+2+nappend]=omega*omega; staFunc->NormalizeVector(L); } polufinal: if (Ret) { v1=staFunc->DestroyVector(v1); v2=staFunc->DestroyVector(v2); D=staFunc->DestroyMatrix(D); DE=staFunc->DestroyMatrix(DE); Term(FALSE); return 1; } } } if (JacRHS(x0,C1)) return 1; /* for safety */ staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->ClearMatrix(BP); staFunc->MultiplyMatrixMatrix(A,A,BP); for (i=0;i= 2) { staFunc->ClearMatrix(BP1); staFunc->CopyMatrixBlock(BP,0,0,BP1,0,0,nphase,nphase); for (i=0; i t) { p = i; q = j; t = fabs(BP1[i][j]); } } } if (k!=q) { /* switch colom k with colom q */ for (i=0;i<=nphase;i++) { t=BP1[i][k]; BP1[i][k]=BP1[i][q]; BP1[i][q]=t; } } if (k!=p) { /* switch row k with row p */ for (i=0;i<=nphase;i++) { t=BP1[k][i]; BP1[k][i]=BP1[p][i]; BP1[p][i]=t; } } if (fabs(BP1[k][k]) == 0.0) { staUtil->ErrMessage("(A^2 + k I) has rank defect at least three"); Term(FALSE); return 1; } kp1=k+1; for (i=kp1;iClearVector(V1); staFunc->ClearVector(V2); staFunc->ClearVector(W1); staFunc->ClearVector(W2); staFunc->CopyMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(VB1); VB1[nphase]=1.0; staFunc->Solve(BT,VB1); staFunc->CopyVector(VB1,V1); staFunc->ClearVector(VB2); VB2[nphase+1]=1.0; staFunc->Solve(BT,VB2); staFunc->CopyVector(VB2,V2); /* solve transposed systems to obtain new values for W1 and W2*/ staFunc->TransposeMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(PP1); PP1[nphase]=1.0; staFunc->Solve(BT,PP1); staFunc->CopyVector(PP1,W1); staFunc->ClearVector(PP2); PP2[nphase+1]=1.0; staFunc->Solve(BT,PP2); staFunc->CopyVector(PP2,W2); for (i=0;iClearVector(V1); staFunc->ClearVector(V2); staFunc->ClearVector(W1); staFunc->ClearVector(W2); V1[0]=1.0; W1[0]=1.0; W2[nphase+1]=1.0; } staFunc->CopyMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(VB1); VB1[nphase]=1.0; staFunc->Solve(BT,VB1); staFunc->ClearVector(VB2); VB2[nphase+1]=1.0; staFunc->Solve(BT,VB2); staFunc->TransposeMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(PP1); PP1[nphase]=1.0; staFunc->Solve(BT,PP1); staFunc->ClearVector(PP2); PP2[nphase+1]=1.0; staFunc->Solve(BT,PP2); /* compute derivatives of g11,g12,g21,g22 */ staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); if (stagenData->Der2) { CopyToX(x0); CopyToP(x0); H=stagenData->Der2(X,P,&zero); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->MultiplyMatrixVector(A,HV1,HV1); staFunc->MultiplyMatrixVector(A,HV2,HV2); staFunc->TransposeMatrix(A,B); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); staFunc->MultiplyMatrixVector(B,HW1,HW1); staFunc->MultiplyMatrixVector(B,HW2,HW2); S=staFunc->CreateMatrix(2,2); for (l=0; lDestroyMatrix(S); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); G0[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV1); G1[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV2); G2[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV1); G3[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV2); } else { staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->MultiplyMatrixVector(A,HV1,HV1); staFunc->MultiplyMatrixVector(A,HV2,HV2); staFunc->TransposeMatrix(A,B); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); staFunc->MultiplyMatrixVector(B,HW1,HW1); staFunc->MultiplyMatrixVector(B,HW2,HW2); Jac1=staFunc->CreateMatrix(nphase,nphase); Jac2=staFunc->CreateMatrix(nphase,nphase); for (l=0; lClearVector(xp); staFunc->CopyVectorElements(x0,0,xp,0,nphase+2); xp[l]+=_dx; staFunc->ClearVector(xp1); staFunc->CopyVectorElements(x0,0,xp1,0,nphase+2); xp1[l]-=_dx; JacRHS(xp,C1); staFunc->CopyMatrixBlock(C1,0,0,Jac1,0,0,nphase,nphase); JacRHS(xp1,C1); staFunc->CopyMatrixBlock(C1,0,0,Jac2,0,0,nphase,nphase); /* g11 */ staFunc->CopyVectorElements(VB1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G0[l]=-staFunc->ScalProd(HW1,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G0[l]-=staFunc->ScalProd(K3,HV1); /* g12 */ staFunc->CopyVectorElements(VB2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G1[l]=-staFunc->ScalProd(HW1,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G1[l]-=staFunc->ScalProd(K3,HV2); /* g21 */ staFunc->CopyVectorElements(VB1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G2[l]=-staFunc->ScalProd(HW2,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G2[l]-=staFunc->ScalProd(K3,HV1); /* g22 */ staFunc->CopyVectorElements(VB2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G3[l]=-staFunc->ScalProd(HW2,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G3[l]-=staFunc->ScalProd(K3,HV2); } Jac1=staFunc->DestroyMatrix(Jac1); Jac2=staFunc->DestroyMatrix(Jac2); G0[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV1); G1[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV2); G2[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV1); G3[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV2); } staFunc->ClearMatrix(BT); staFunc->ClearMatrix(JM1); staFunc->CopyMatrixBlock(C1,0,0,BT,0,0,nphase,nphase+2); /* placement of ones determined by partial pivoting */ PIVOT=staFunc->CreateVector(nphase+2); staFunc->ClearVector(PIVOT); staFunc->TransposeMatrix(BT,JM1); PIVOT[nphase+1]=1.0; for (k=0;k<=nphase;k++) { kp1=k+1; p = k; for (i=kp1;i<=nphase+1;i++) { if (fabs(JM1[i][k]) > fabs(JM1[p][k])) p = i; } PIVOT[k]=p; if (p!=k) PIVOT[nphase+1]=-PIVOT[nphase+1]; t=JM1[p][k]; JM1[p][k]=JM1[k][k]; JM1[k][k]=t; if (t==0.0) continue; for (i=kp1;i<=nphase+1;i++) JM1[i][k]=-JM1[i][k]/t; for (j=kp1;j<=nphase+1;j++) { t=JM1[p][j]; JM1[p][j]=JM1[k][j]; JM1[k][j]=t; if (t!=0.0) for (i=kp1;i<=nphase+1;i++) JM1[i][j]+=JM1[i][k]*t; } } staFunc->ClearVector(TV1); for (i=0;iDestroyVector(PIVOT); /* end of pivot information */ staFunc->ClearMatrix(JM1); staFunc->CopyMatrixBlock(C1,0,0,JM1,0,0,nphase,nphase+2); JM1[nphase][p]=1.0; JM1[nphase+1][q]=1.0; staFunc->ClearVector(TV1); TV1[nphase]=1.0; Ret=staFunc->Decomp(JM1); if (Ret) return 1; staFunc->Solve(JM1,TV1); staFunc->ClearMatrix(TVM); for (i=0;iClearVector(TV1); TV1[nphase+1]=1.0; staFunc->Solve(JM1,TV1); for (i=0;iTransposeMatrix(TVM,TVMT); staFunc->MultiplyMatrixMatrix(TVMT,TVM,HM); Ret=staFunc->Decomp(HM); if (Ret) return 1; staFunc->MultiplyMatrixVector(TVMT,G0,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G0); staFunc->MultiplyMatrixVector(TVMT,G1,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G1); staFunc->MultiplyMatrixVector(TVMT,G2,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G2); staFunc->MultiplyMatrixVector(TVMT,G3,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G3); /* Determine component with maximum norm => gives m1,n1 */ G[0][0] = staFunc->NormOfVector(G0); G[0][1] = staFunc->NormOfVector(G1); G[1][0] = staFunc->NormOfVector(G2); G[1][1] = staFunc->NormOfVector(G3); m1=0; n1=0; if (G[0][1] > G[m1][n1]) n1=1; if (G[1][0] > G[m1][n1]) { m1=1; n1=0; } if (G[1][1] > G[m1][n1]) { m1=1; n1=1; } for (i=0;iNormalizeVector(G0); /* G1,G2,G3 are projected on G0 */ for (j=0;j<4;j++) { if (j != (2*m1+n1)) { for (i=0;iScalProd(G1,G0); staFunc->AddScaledVector(G1,G0,-scal,G1); for (i=0;i gives m2,n2 */ for (i=0;iNormOfVector(G0); G[0][1] = staFunc->NormOfVector(G1); G[1][0] = staFunc->NormOfVector(G2); G[1][1] = staFunc->NormOfVector(G3); G[m1][n1]=-1.0; m2=0; n2=0; if (G[0][1] > G[m2][n2]) n2=1; if (G[1][0] > G[m2][n2]) { m2=1; n2=0; } if (G[1][1] > G[m2][n2]) { m2=1; n2=1; } /* Compute a tangent vector to the Hopf curve */ h2_DefHopf(x0,v0); Ret=h2_JacDefHopf(x0,D); for (i=0; iAppendMatrixVector(D,v1,DE); Ret=staFunc->Decomp(DE); if (Ret == 0) { v2[nphase+2+nappend]=stagenData->Forward ? 1 : -1; staFunc->Solve(DE,v2); staFunc->CopyVector(v2,v0); break; }; v1[i]=0.0; }; if (Ret) { staUtil->ErrMessage ("Unable to find a tangent vector to the Hopf curve"); v1=staFunc->DestroyVector(v1); v2=staFunc->DestroyVector(v2); D=staFunc->DestroyMatrix(D); DE=staFunc->DestroyMatrix(DE); Term(FALSE); return 1; } staFunc->NormalizeVector(v0); v1=staFunc->DestroyVector(v1); v2=staFunc->DestroyVector(v2); D=staFunc->DestroyMatrix(D); DE=staFunc->DestroyMatrix(DE); } /* Initialize continuation */ ((ContDataPtr)stagenptr->GenPtr)->X0=x0; ((ContDataPtr)stagenptr->GenPtr)->V0=v0; ((ContDataPtr)stagenptr->GenPtr)->ActiveSing=staData->ActiveSing; ((ContDataPtr)stagenptr->GenPtr)->ActiveUserTest=staData->ActiveUserTest; ((ContDataPtr)stagenptr->GenPtr)->TypeTestFuns=staData->TypeTestFunc; ((ContDataPtr)stagenptr->GenPtr)->ZeroTestFuns=staData->ZeroTestFunc; /* That's all */ return 0; } /*************************************************************/ /* Defining conditions for the eigenvector: (A*A+kappa*E)V=0 */ Entry(Int2) DefEig (Vector x, Vector y) { if (JacRHS(x,C1)) return 1; /* A,C1,V,W - global */ staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->CopyVectorElements(x, nphase+2+nappend, V, 0, nphase); staFunc->MultiplyMatrixVector(A,V,W); staFunc->MultiplyMatrixVector(A,W,V); staFunc->CopyVectorElements(x, nphase+2+nappend, W, 0, nphase); staFunc->AddScaledVector(V,W,x[2*nphase+2+nappend],V); staFunc->CopyVectorElements(V, 0, y, 0, nphase); return 0; } /*******************/ /* dh -> h2 Starter */ Entry(Int2) dh_h2_Starter(StagenDataPtr stagenptr) { Int2 i,j,k,kp1,l,Ret,p,q; Vector v1,v2; Matrix D,DE,S; VAR t,scal; if (!stagenptr) { /* this is request to terminate */ Term(TRUE); return 0; } if (Init(stagenptr)) return 1; /* initialize local variables */ /* Create the first point for generator */ if (stagenptr->ContDataGenLen) { /* "Extension" */ MemCopy(V1,stagenptr->ContDataStaPtr,sizeof(FloatHi)*(nphase+2)); MemCopy(W1,stagenptr->ContDataStaPtr+(nphase+2),sizeof(FloatHi)*(nphase+2)); MemCopy(V2,stagenptr->ContDataStaPtr+2*(nphase+2),sizeof(FloatHi)*(nphase+2)); MemCopy(W2,stagenptr->ContDataStaPtr+3*(nphase+2),sizeof(FloatHi)*(nphase+2)); MemFree(stagenptr->ContDataStaPtr); stagenptr->ContDataStaLen=0; MemCopy(x0,staData->X,sizeof(FloatHi)*nphase); for (i=0; i<2+nappend; i++) x0[nphase+i]=staData->P[active[i]]; x0[nphase+2+nappend]=BifVector[2*nphase]; } else { v1=staFunc->CreateVector(nphase+3+nappend); v2=staFunc->CreateVector(nphase+3+nappend); D=staFunc->CreateMatrix(nphase+2+nappend,nphase+3+nappend); DE=staFunc->CreateMatrix(nphase+3+nappend,nphase+3+nappend); MemCopy(x0,staData->X,sizeof(FloatHi)*nphase); for (i=0; i<2+nappend; i++) x0[nphase+i]=staData->P[active[i]]; JacRHS(x0,C1); if (stagenData->BifDataInOk) { /* Message(MSG_OK,"Eigenvectors exist"); */ /* bifurcation data exist */ /* BifVector[0] to [nphase-1]: Hopf continuation vector */ /* BifVector[nphase] to [2*nphase-1]: global vector L */ /* BifVector[2*nphase]: -lambda*lambda */ MemCopy(V,stagenData->BifDataInPtr,sizeof(FloatHi)*nphase); MemCopy(L,stagenData->BifDataInPtr+nphase,sizeof(FloatHi)*nphase); x0[nphase+2+nappend]=stagenData->BifDataInPtr[2*nphase]; } else { /* bifurcation data to be computed */ /* Message(MSG_OK,"Eigenvectors have to be computed"); */ Ret=0; {int i,j,imin,jmin; double r,s,d,phi,beta,gamma,lambda,omega; staFunc->CopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); staFunc->EigenValues(AA,Re,Im); /* find the pair of eigenvalues with zero sum */ s=fabs(Re[0]+Re[1])+fabs(Im[0]+Im[1]); imin=0; jmin=1; for (i=0; iCopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); for (i=0; iSngSolve(AA,_eps,V) != nphase-1) { staUtil->ErrMessage("Unable to find (+)eigenvector at the initial point"); Ret=1; goto polufinal; } staFunc->CopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); for (i=0; iSngSolve(AA,_eps,W) != nphase-1) { staUtil->ErrMessage("Unable to find (-)eigenvector at the initial point"); Ret=1; goto polufinal; } staFunc->NormalizeVector(V); staFunc->NormalizeVector(W); staFunc->AddScaledVector(V,W,-1.0,L); staFunc->NormalizeVector(L); staFunc->AddScaledVector(V,W,1.0,V); staFunc->NormalizeVector(V); x0[nphase+2+nappend]=-lambda*lambda; } else { /* Hopf */ lambda=Re[imin]; omega=fabs(Im[imin]); /* compute real and imaginary part of the eigenvector */ staFunc->CopyMatrixBlock(C1,0,0,BB,0,0,nphase,nphase); staFunc->CopyMatrixBlock(C1,0,0,BB,nphase,nphase,nphase,nphase); for (i=0; iSngSolve(BB,_eps,VV) != 2*nphase-2) { staUtil->ErrMessage("Unable to find eigenvectors at the initial point"); Ret=1; goto polufinal; } staFunc->CopyVectorElements(VV,0,V,0,nphase); staFunc->CopyVectorElements(VV,nphase,W,0,nphase); /* normalize real and imaginary parts */ d=staFunc->ScalProd(V,V); s=staFunc->ScalProd(W,W); r=staFunc->ScalProd(V,W); if (r!=0) { beta=sqrt(fabs(4.0*r*r+(s-d)*(s-d))); if (2.0*fabs(2.0*r) > beta) { phi=PID2-atan(fabs((s-d)/(2.0*r))); } else { phi=atan(fabs(2.0*r/(s-d))); } if (r < 0) phi=PI-phi; if ((s-d) < 0) phi=-phi; phi/=2.0; for (i=0; iCopyVectorElements(VV,0,V,0,nphase); d=staFunc->ScalProd(V,V); gamma=1/sqrt(d); staFunc->MultiplyVectorScalar(VV,gamma,VV); staFunc->CopyVectorElements(VV,0,V,0,nphase); staFunc->CopyVectorElements(VV,nphase,L,0,nphase); x0[nphase+2+nappend]=omega*omega; staFunc->NormalizeVector(L); } polufinal: if (Ret) { v1=staFunc->DestroyVector(v1); v2=staFunc->DestroyVector(v2); D=staFunc->DestroyMatrix(D); DE=staFunc->DestroyMatrix(DE); Term(FALSE); return 1; } } } if (JacRHS(x0,C1)) return 1; /* for safety */ staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->ClearMatrix(BP); staFunc->MultiplyMatrixMatrix(A,A,BP); for (i=0;i= 2) { staFunc->ClearMatrix(BP1); staFunc->CopyMatrixBlock(BP,0,0,BP1,0,0,nphase,nphase); for (i=0; i t) { p = i; q = j; t = fabs(BP1[i][j]); } } } if (k!=q) { /* switch colom k with colom q */ for (i=0;i<=nphase;i++) { t=BP1[i][k]; BP1[i][k]=BP1[i][q]; BP1[i][q]=t; } } if (k!=p) { /* switch row k with row p */ for (i=0;i<=nphase;i++) { t=BP1[k][i]; BP1[k][i]=BP1[p][i]; BP1[p][i]=t; } } if (fabs(BP1[k][k]) == 0.0) { staUtil->ErrMessage("(A^2 + k I) has rank defect at least three"); Term(FALSE); return 1; } kp1=k+1; for (i=kp1;iClearVector(V1); staFunc->ClearVector(V2); staFunc->ClearVector(W1); staFunc->ClearVector(W2); staFunc->CopyMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(VB1); VB1[nphase]=1.0; staFunc->Solve(BT,VB1); staFunc->CopyVector(VB1,V1); staFunc->ClearVector(VB2); VB2[nphase+1]=1.0; staFunc->Solve(BT,VB2); staFunc->CopyVector(VB2,V2); /* solve transposed systems to obtain new values for W1 and W2*/ staFunc->TransposeMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(PP1); PP1[nphase]=1.0; staFunc->Solve(BT,PP1); staFunc->CopyVector(PP1,W1); staFunc->ClearVector(PP2); PP2[nphase+1]=1.0; staFunc->Solve(BT,PP2); staFunc->CopyVector(PP2,W2); for (i=0;iClearVector(V1); staFunc->ClearVector(V2); staFunc->ClearVector(W1); staFunc->ClearVector(W2); V1[0]=1.0; W1[0]=1.0; W2[nphase+1]=1.0; } staFunc->CopyMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(VB1); VB1[nphase]=1.0; staFunc->Solve(BT,VB1); staFunc->ClearVector(VB2); VB2[nphase+1]=1.0; staFunc->Solve(BT,VB2); staFunc->TransposeMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(PP1); PP1[nphase]=1.0; staFunc->Solve(BT,PP1); staFunc->ClearVector(PP2); PP2[nphase+1]=1.0; staFunc->Solve(BT,PP2); /* compute derivatives of g11,g12,g21,g22 */ staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); if (stagenData->Der2) { CopyToX(x0); CopyToP(x0); H=stagenData->Der2(X,P,&zero); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->MultiplyMatrixVector(A,HV1,HV1); staFunc->MultiplyMatrixVector(A,HV2,HV2); staFunc->TransposeMatrix(A,B); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); staFunc->MultiplyMatrixVector(B,HW1,HW1); staFunc->MultiplyMatrixVector(B,HW2,HW2); S=staFunc->CreateMatrix(2,2); for (l=0; lDestroyMatrix(S); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); G0[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV1); G1[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV2); G2[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV1); G3[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV2); } else { staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->MultiplyMatrixVector(A,HV1,HV1); staFunc->MultiplyMatrixVector(A,HV2,HV2); staFunc->TransposeMatrix(A,B); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); staFunc->MultiplyMatrixVector(B,HW1,HW1); staFunc->MultiplyMatrixVector(B,HW2,HW2); Jac1=staFunc->CreateMatrix(nphase,nphase); Jac2=staFunc->CreateMatrix(nphase,nphase); for (l=0; lClearVector(xp); staFunc->CopyVectorElements(x0,0,xp,0,nphase+2); xp[l]+=_dx; staFunc->ClearVector(xp1); staFunc->CopyVectorElements(x0,0,xp1,0,nphase+2); xp1[l]-=_dx; JacRHS(xp,C1); staFunc->CopyMatrixBlock(C1,0,0,Jac1,0,0,nphase,nphase); JacRHS(xp1,C1); staFunc->CopyMatrixBlock(C1,0,0,Jac2,0,0,nphase,nphase); /* g11 */ staFunc->CopyVectorElements(VB1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G0[l]=-staFunc->ScalProd(HW1,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G0[l]-=staFunc->ScalProd(K3,HV1); /* g12 */ staFunc->CopyVectorElements(VB2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G1[l]=-staFunc->ScalProd(HW1,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G1[l]-=staFunc->ScalProd(K3,HV2); /* g21 */ staFunc->CopyVectorElements(VB1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G2[l]=-staFunc->ScalProd(HW2,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G2[l]-=staFunc->ScalProd(K3,HV1); /* g22 */ staFunc->CopyVectorElements(VB2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G3[l]=-staFunc->ScalProd(HW2,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G3[l]-=staFunc->ScalProd(K3,HV2); } Jac1=staFunc->DestroyMatrix(Jac1); Jac2=staFunc->DestroyMatrix(Jac2); G0[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV1); G1[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV2); G2[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV1); G3[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV2); } staFunc->ClearMatrix(BT); staFunc->ClearMatrix(JM1); staFunc->CopyMatrixBlock(C1,0,0,BT,0,0,nphase,nphase+2); /* placement of ones determined by partial pivoting */ PIVOT=staFunc->CreateVector(nphase+2); staFunc->ClearVector(PIVOT); staFunc->TransposeMatrix(BT,JM1); PIVOT[nphase+1]=1.0; for (k=0;k<=nphase;k++) { kp1=k+1; p = k; for (i=kp1;i<=nphase+1;i++) { if (fabs(JM1[i][k]) > fabs(JM1[p][k])) p = i; } PIVOT[k]=p; if (p!=k) PIVOT[nphase+1]=-PIVOT[nphase+1]; t=JM1[p][k]; JM1[p][k]=JM1[k][k]; JM1[k][k]=t; if (t==0.0) continue; for (i=kp1;i<=nphase+1;i++) JM1[i][k]=-JM1[i][k]/t; for (j=kp1;j<=nphase+1;j++) { t=JM1[p][j]; JM1[p][j]=JM1[k][j]; JM1[k][j]=t; if (t!=0.0) for (i=kp1;i<=nphase+1;i++) JM1[i][j]+=JM1[i][k]*t; } } staFunc->ClearVector(TV1); for (i=0;iDestroyVector(PIVOT); /* end of pivot information */ staFunc->ClearMatrix(JM1); staFunc->CopyMatrixBlock(C1,0,0,JM1,0,0,nphase,nphase+2); JM1[nphase][p]=1.0; JM1[nphase+1][q]=1.0; staFunc->ClearVector(TV1); TV1[nphase]=1.0; Ret=staFunc->Decomp(JM1); if (Ret) return 1; staFunc->Solve(JM1,TV1); staFunc->ClearMatrix(TVM); for (i=0;iClearVector(TV1); TV1[nphase+1]=1.0; staFunc->Solve(JM1,TV1); for (i=0;iTransposeMatrix(TVM,TVMT); staFunc->MultiplyMatrixMatrix(TVMT,TVM,HM); Ret=staFunc->Decomp(HM); if (Ret) return 1; staFunc->MultiplyMatrixVector(TVMT,G0,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G0); staFunc->MultiplyMatrixVector(TVMT,G1,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G1); staFunc->MultiplyMatrixVector(TVMT,G2,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G2); staFunc->MultiplyMatrixVector(TVMT,G3,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G3); /* Determine component with maximum norm => gives m1,n1 */ G[0][0] = staFunc->NormOfVector(G0); G[0][1] = staFunc->NormOfVector(G1); G[1][0] = staFunc->NormOfVector(G2); G[1][1] = staFunc->NormOfVector(G3); m1=0; n1=0; if (G[0][1] > G[m1][n1]) n1=1; if (G[1][0] > G[m1][n1]) { m1=1; n1=0; } if (G[1][1] > G[m1][n1]) { m1=1; n1=1; } for (i=0;iNormalizeVector(G0); /* G1,G2,G3 are projected on G0 */ for (j=0;j<4;j++) { if (j != (2*m1+n1)) { for (i=0;iScalProd(G1,G0); staFunc->AddScaledVector(G1,G0,-scal,G1); for (i=0;i gives m2,n2 */ for (i=0;iNormOfVector(G0); G[0][1] = staFunc->NormOfVector(G1); G[1][0] = staFunc->NormOfVector(G2); G[1][1] = staFunc->NormOfVector(G3); G[m1][n1]=-1.0; m2=0; n2=0; if (G[0][1] > G[m2][n2]) n2=1; if (G[1][0] > G[m2][n2]) { m2=1; n2=0; } if (G[1][1] > G[m2][n2]) { m2=1; n2=1; } /* Compute a tangent vector to the Hopf curve */ h2_DefHopf(x0,v0); Ret=h2_JacDefHopf(x0,D); for (i=0; iAppendMatrixVector(D,v1,DE); Ret=staFunc->Decomp(DE); if (Ret == 0) { v2[nphase+2+nappend]=stagenData->Forward ? 1 : -1; staFunc->Solve(DE,v2); staFunc->CopyVector(v2,v0); break; }; v1[i]=0.0; }; if (Ret) { staUtil->ErrMessage ("Unable to find a tangent vector to the Hopf curve"); v1=staFunc->DestroyVector(v1); v2=staFunc->DestroyVector(v2); D=staFunc->DestroyMatrix(D); DE=staFunc->DestroyMatrix(DE); Term(FALSE); return 1; } staFunc->NormalizeVector(v0); v1=staFunc->DestroyVector(v1); v2=staFunc->DestroyVector(v2); D=staFunc->DestroyMatrix(D); DE=staFunc->DestroyMatrix(DE); } /* Initialize continuation */ ((ContDataPtr)stagenptr->GenPtr)->X0=x0; ((ContDataPtr)stagenptr->GenPtr)->V0=v0; ((ContDataPtr)stagenptr->GenPtr)->ActiveSing=staData->ActiveSing; ((ContDataPtr)stagenptr->GenPtr)->ActiveUserTest=staData->ActiveUserTest; ((ContDataPtr)stagenptr->GenPtr)->TypeTestFuns=staData->TypeTestFunc; ((ContDataPtr)stagenptr->GenPtr)->ZeroTestFuns=staData->ZeroTestFunc; /* That's all */ return 0; } /************************/ /* Some common routines */ Local(void) CopyToX(Vector vp) { MemCopy(X,vp,sizeof(FloatHi)*nphase); } Local(void) CopyToP(Vector vp) { Int2 i; for (i=0; i<2+nappend; i++) P[active[i]]=vp[nphase+i]; } Local(void) CopyToE(Vector Re, Vector Im) { MemCopy(E,Re,sizeof(FloatHi)*nphase); MemCopy(E+nphase,Im,sizeof(FloatHi)*nphase); } /***********/ /* Decoder */ Entry(Int2) h2_Decoder(VoidPtr vpoint, FloatHiPtr PNTR indirect, Int2 oper) { switch (oper) { case DECODER_INIT: /* point is a pointer to StagenData structure */ case DECODER_TERM: break; case DECODER_DIM: /* returns number of M-points in given G-point */ return 1; default: /* extract next M-point from G-point pointed by vpoint */ /* vpoint is a Vector */ indirect[CL_TIME]=&zero; indirect[CL_PHASE]=X; CopyToX(vpoint); indirect[CL_PAR]=P; CopyToP(vpoint); indirect[CL_EIGV]=E; CopyToE(Re,Im); indirect[CL_TEST]=test; } return 0; } /***************************/ /* User-defined functions */ Entry(Int2) h2_DefUser (Vector x, Vector y) { Int2 i,j; for (i=j=0; iudFunc)(x,i,&y[j]); j++; } return 0; } /*******/ /* RHS */ Entry(Int2) h2_DefRHS (Vector x, Vector y) { CopyToX(x); CopyToP(x); stagenData->Rhs(X,P,&zero,y); return 0; } /*************************************/ /* Defining function for equilibrium */ Local(Int2) DefEquilib (Vector x, Vector y) { h2_DefRHS(x,y); if (nappend) h2_DefUser(x,y+nphase); return 0; } /********************/ /* Jacobian of RHS */ Local(Int2) JacRHS (Vector x, Matrix Jac) { staFunc->CopyVectorElements(x,0,xp,0,nphase+2+nappend); if (stagenData->Der1) { SymJac(xp,Jac); } else { if (staFunc->Der(h2_DefRHS, nphase, xp, dx, Jac)) { return 1; } } return 0; } /***************************************************/ /* Jacobian of defining conditions for equilibria */ Local(Int2) JacDefEquilib (Vector x, Matrix Jac) { staFunc->ClearMatrix(Jac); JacRHS(x, Jac); if (nappend) staFunc->Der(h2_DefUser, nappend, x, dx, Jac+nphase); return 0; } /***************************/ /* Hopf defining condition */ Entry(Int2) h2_DefBord (Vector x, Vector y) { Int2 i,Ret; if (JacRHS(x,C1)) return 1; /* A,C1,BP,V1,V2,W1,W2 - global */ staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->MultiplyMatrixMatrix(A,A,BP); for (i=0;iCopyMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(VB1); VB1[nphase]=1.0; staFunc->Solve(BT,VB1); staFunc->CopyVector(VB1,Q1); staFunc->ClearVector(VB2); VB2[nphase+1]=1.0; staFunc->Solve(BT,VB2); staFunc->CopyVector(VB2,Q2); G[0][0]=VB1[nphase]; /* g11 */ G[1][0]=VB1[nphase+1]; /* g21 */ G[0][1]=VB2[nphase]; /* g12 */ G[1][1]=VB2[nphase+1]; /* g22 */ y[0]=G[m1][n1]; y[1]=G[m2][n2]; return 0; } /***************************************/ /* Jacobian of Hopf defining condition */ Local(Int2) JacDefBord (Vector x, Matrix Jac) { Int2 i,j,l,Ret; Matrix S; staFunc->ClearMatrix(Jac); if (stagenData->Der2) { /* creation of DefBord Jacobian matrix Jac using symbolic second derivatives */ CopyToX(x); CopyToP(x); H=stagenData->Der2(X,P,&zero); /* take global vectors VB1 and VB2 computed in DefBord */ /* compute PP1 cand PP2 using the tranposed squared Jacobian */ staFunc->TransposeMatrix(BP,BT); Ret=staFunc->Decomp(BT); if (Ret) return 1; staFunc->ClearVector(PP1); PP1[nphase]=1.0; staFunc->Solve(BT,PP1); staFunc->ClearVector(PP2); PP2[nphase+1]=1.0; staFunc->Solve(BT,PP2); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,VB,0,nphase); staFunc->MultiplyMatrixVector(A,VB,HV1); staFunc->CopyVectorElements(VB2,0,VB,0,nphase); staFunc->MultiplyMatrixVector(A,VB,HV2); staFunc->TransposeMatrix(A,B); staFunc->CopyVectorElements(PP1,0,VB,0,nphase); staFunc->MultiplyMatrixVector(B,VB,HW1); staFunc->CopyVectorElements(PP2,0,VB,0,nphase); staFunc->MultiplyMatrixVector(B,VB,HW2); S=staFunc->CreateMatrix(2,2); for (l=0; lDestroyMatrix(S); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); G0[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV1); G1[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV2); G2[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV1); G3[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV2); GM1[nphase+2+nappend][0]=G0[nphase+2+nappend]; GM1[nphase+2+nappend][1]=G1[nphase+2+nappend]; GM1[nphase+2+nappend][2]=G2[nphase+2+nappend]; GM1[nphase+2+nappend][3]=G3[nphase+2+nappend]; Jac[0][nphase+2+nappend]=GM1[nphase+2+nappend][2*m1+n1]; Jac[1][nphase+2+nappend]=GM1[nphase+2+nappend][2*m2+n2]; } else { if (staFunc->Der(h2_DefBord, 2, x, dx, Jac)) return 1; } return(0); } /*************************************/ /* Defining functions for Hopf curve */ Entry(Int2) h2_DefHopf (Vector x, Vector y) { if (DefEquilib(x,y)) return 1; if (h2_DefBord(x,y+nphase+nappend)) return 2; return 0; } /*************************************************/ /* Jacobian of defining functions for Hopf curve */ Entry(Int2) h2_JacDefHopf (Vector x, Matrix Jac) { staFunc->ClearMatrix(Jac); if (JacDefEquilib(x,C)) return 1; staFunc->CopyMatrixBlock(C,0,0,Jac,0,0,nphase+nappend,nphase+2+nappend); if (JacDefBord(x,C2)) return 2; staFunc->CopyMatrixBlock(C2,0,0,Jac,nphase+nappend,0,2,nphase+3+nappend); return 0; } /*********************/ /* Default processor */ Entry(Int2) h2_ProcDefault(Int2 Ind, Vector x, Vector v1, CharPtr *msg) { /* Assert: v1==NULL, msg==NULL, Ind==0 */ h2_DefHopf (x, xp1); /* compute global Q for test functions */ MemCopy(x1,x,sizeof(FloatHi)*(nphase+3+nappend)); /* for the last point */ if (!staData->eigcomp) return 0; if (JacRHS(x,C1)) return (1); staFunc->EigenValues(C1,Re,Im); return 0; } /*********************************/ /* Test and processing functions */ typedef EntryPtr(Int2,TestFuncPtr,(Vector xx, Vector vv, FloatHiPtr TV)); typedef EntryPtr(Int2,ProcFuncPtr,(Int2 Ind, Vector xx, Vector vv, CharPtr *msg)); Local(Char) buf[80]; /* first Lyapunov coefficient */ Entry(Int2) h2_TestBautin(Vector x, Vector v, FloatHiPtr res) { VAR beta,gamma,phi,omega,s,r,h,d; Int2 i,j,k,l,Ret; /* reevaluate Jacobian */ if (JacRHS (x,C1)) { return 1; } staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); if (x[nphase+2+nappend] <= 0.0) { /* neutral saddle */ *res=777; } else { /* Hopf */ omega=sqrt(x[nphase+2+nappend]); staFunc->CopyVectorElements(VB1,0,V,0,nphase); staFunc->CopyVectorElements(VB2,0,W,0,nphase); staFunc->MultiplyMatrixVector(A,V,W); staFunc->MultiplyVectorScalar(W,-1.0/omega,W); /* normalize real and imaginary parts =1, =0 */ d=staFunc->ScalProd(V,V); s=staFunc->ScalProd(W,W); r=staFunc->ScalProd(V,W); if (r==0) { staFunc->CopyVectorElements(V,0,VV,0,nphase); staFunc->CopyVectorElements(W,0,VV,nphase,nphase); } else { beta=sqrt(fabs(4.0*r*r+(s-d)*(s-d))); if (2.0*fabs(2.0*r) > beta) { phi=PID2-atan(fabs((s-d)/(2.0*r))); } else { phi=atan(fabs(2.0*r/(s-d))); } if (r < 0) phi=PI-phi; if ((s-d) < 0) phi=-phi; phi/=2.0; for (i=0; iCopyVectorElements(VV,0,W,0,nphase); d=staFunc->ScalProd(W,W); gamma=1/sqrt(d); staFunc->MultiplyVectorScalar(VV,gamma,VV); staFunc->CopyVectorElements(VV,0,Qr,0,nphase); staFunc->CopyVectorElements(VV,nphase,Qi,0,nphase); /* compute real and imaginary part of the adjoint eigenvector */ staFunc->CopyMatrixBlock(A,0,0,AA,0,0,nphase,nphase); staFunc->TransposeMatrix(AA,AT); staFunc->ClearMatrix(BB); staFunc->CopyMatrixBlock(AT,0,0,BB,0,0,nphase,nphase); staFunc->CopyMatrixBlock(AT,0,0,BB,nphase,nphase,nphase,nphase); for (i=0; iSngSolve(BB,_eps,VV) != 2*nphase-2) { *res=1.0; test[0]=*res; return 1; } staFunc->CopyVectorElements(VV,0,Pr,0,nphase); staFunc->CopyVectorElements(VV,nphase,Pi,0,nphase); /* normalize the adjoint vector +=1, -=0 */ d=staFunc->ScalProd(Pr,Qr); s=staFunc->ScalProd(Pi,Qi); r=staFunc->ScalProd(Pr,Qi); h=staFunc->ScalProd(Pi,Qr); /* phi=atan((r-h)/(d+s)); */ beta=sqrt(fabs((r-h)*(r-h)+(d+s)*(d+s))); if (2.0*fabs(d+s) > beta) { phi=atan(fabs((r-h)/(d+s))); } else { phi=PID2-atan(fabs((d+s)/(r-h))); } if ((r-h) < 0) phi=PI-phi; if ((d+s) < 0) phi=-phi; for (i=0; iCopyVectorElements(VV,0,Pr,0,nphase); staFunc->CopyVectorElements(VV,nphase,Pi,0,nphase); r=staFunc->ScalProd(Pr,Qr); h=staFunc->ScalProd(Pi,Qi); beta=1.0/(r+h); staFunc->MultiplyVectorScalar(Pr,beta,Pr); staFunc->MultiplyVectorScalar(Pi,beta,Pi); /* compute the first Lyapunov coefficient L1 */ /* quadratic terms */ if (stagenData->Der2) { CopyToX(x); CopyToP(x); H=stagenData->Der2(X,P,&zero); /* The 1st elem must be time! */ for (i=0; iClearVector(v0); staFunc->CopyVectorElements(Qr,0,v0,0,nphase); h2_DefRHS(x,aa); staFunc->AddScaledVector(x,v0,ddx,x1); h2_DefRHS(x1,bb); staFunc->AddScaledVector(x,v0,-ddx,x1); h2_DefRHS(x1,cc); staFunc->AddScaledVector(bb,aa,-2.0,rr); staFunc->AddScaledVector(rr,cc,1.0,aa); staFunc->MultiplyVectorScalar(aa,1.0/ddx/ddx,aa); staFunc->CopyVectorElements(Qi,0,v0,0,nphase); h2_DefRHS(x,rr); staFunc->AddScaledVector(x,v0,ddx,x1); h2_DefRHS(x1,bb); staFunc->AddScaledVector(x,v0,-ddx,x1); h2_DefRHS(x1,cc); staFunc->AddScaledVector(bb,rr,-2.0,rr); staFunc->AddScaledVector(rr,cc,1.0,bb); staFunc->MultiplyVectorScalar(bb,1.0/ddx/ddx,bb); staFunc->AddScaledVector(Qr,Qi,1.0,rr); staFunc->CopyVectorElements(rr,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,rr); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,cc); staFunc->AddScaledVector(rr,cc,1.0,dd); staFunc->AddScaledVector(Qr,Qi,-1.0,rr); staFunc->CopyVectorElements(rr,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,rr); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,cc); staFunc->AddScaledVector(rr,cc,1.0,rr); staFunc->AddScaledVector(dd,rr,-1.0,cc); staFunc->MultiplyVectorScalar(cc,-2.0/ddx/ddx,cc); /* cc=-2*c */ } staFunc->CopyMatrixBlock(A,0,0,AA,0,0,nphase,nphase); staFunc->AddScaledVector(aa,bb,1.0,rr); Ret=staFunc->Decomp(AA); if (Ret != 0) { /* Zero-Hopf case */ *res=777.0; test[0]=*res; return 0; } staFunc->Solve(AA,rr); staFunc->MultiplyVectorScalar(rr,-2.0,rr); /* rr=-2*r */ staFunc->AddScaledVector(bb,aa,-1.0,dd); staFunc->CopyVectorElements(dd,0,VV,0,nphase); staFunc->CopyVectorElements(cc,0,VV,nphase,nphase); staFunc->ClearMatrix(BB); staFunc->CopyMatrixBlock(A,0,0,BB,0,0,nphase,nphase); staFunc->CopyMatrixBlock(A,0,0,BB,nphase,nphase,nphase,nphase); for (i=0; iDecomp(BB); if (Ret != 0) { *res=1.0; test[0]=*res; return 3; } staFunc->Solve(BB,VV); staFunc->CopyVectorElements(VV,0,aa,0,nphase); /* aa=Sr */ staFunc->CopyVectorElements(VV,nphase,bb,0,nphase); /* bb=Si */ if (stagenData->Der2) { for (s=0.0,i=0; iClearVector(v0); /* -2* */ staFunc->AddScaledVector(Qr,rr,1.0,dd); staFunc->CopyVectorElements(dd,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,dd); s=staFunc->ScalProd(Pr,dd); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,dd); s+=staFunc->ScalProd(Pr,dd); staFunc->AddScaledVector(Qr,rr,-1.0,dd); staFunc->CopyVectorElements(dd,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,dd); s-=staFunc->ScalProd(Pr,dd); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,dd); s-=staFunc->ScalProd(Pr,dd); /* -2* */ staFunc->AddScaledVector(Qi,rr,1.0,dd); staFunc->CopyVectorElements(dd,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,dd); s+=staFunc->ScalProd(Pi,dd); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,dd); s+=staFunc->ScalProd(Pi,dd); staFunc->AddScaledVector(Qi,rr,-1.0,dd); staFunc->CopyVectorElements(dd,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,dd); s-=staFunc->ScalProd(Pi,dd); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,dd); s-=staFunc->ScalProd(Pi,dd); /* */ staFunc->AddScaledVector(Qr,aa,1.0,dd); staFunc->CopyVectorElements(dd,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,dd); s+=staFunc->ScalProd(Pr,dd); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,dd); s+=staFunc->ScalProd(Pr,dd); staFunc->AddScaledVector(Qr,aa,-1.0,dd); staFunc->CopyVectorElements(dd,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,dd); s-=staFunc->ScalProd(Pr,dd); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,dd); s-=staFunc->ScalProd(Pr,dd); /* */ staFunc->AddScaledVector(Qi,bb,1.0,dd); staFunc->CopyVectorElements(dd,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,dd); s+=staFunc->ScalProd(Pr,dd); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,dd); s+=staFunc->ScalProd(Pr,dd); staFunc->AddScaledVector(Qi,bb,-1.0,dd); staFunc->CopyVectorElements(dd,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,dd); s-=staFunc->ScalProd(Pr,dd); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,dd); s-=staFunc->ScalProd(Pr,dd); /* */ staFunc->AddScaledVector(Qr,bb,1.0,dd); staFunc->CopyVectorElements(dd,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,dd); s+=staFunc->ScalProd(Pi,dd); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,dd); s+=staFunc->ScalProd(Pi,dd); staFunc->AddScaledVector(Qr,bb,-1.0,dd); staFunc->CopyVectorElements(dd,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,dd); s-=staFunc->ScalProd(Pi,dd); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,dd); s-=staFunc->ScalProd(Pi,dd); /* */ staFunc->AddScaledVector(Qi,aa,1.0,dd); staFunc->CopyVectorElements(dd,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,dd); s-=staFunc->ScalProd(Pi,dd); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,dd); s-=staFunc->ScalProd(Pi,dd); staFunc->AddScaledVector(Qi,aa,-1.0,dd); staFunc->CopyVectorElements(dd,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,dd); s+=staFunc->ScalProd(Pi,dd); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,dd); s+=staFunc->ScalProd(Pi,dd); s/=(ddx*ddx); } /* cubic terms */ if (stagenData->Der3) { CopyToX(x); CopyToP(x); D3=stagenData->Der3(X,P,&zero); /* The 1st elem must be time! */ for (r=0.0,i=0; iClearVector(v0); staFunc->CopyVectorElements(Qr,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,3*dddx/2,x1); h2_DefRHS(x1,dd); r=(2.0/3.0)*staFunc->ScalProd(Pr,dd); staFunc->AddScaledVector(x,v0,dddx/2,x1); h2_DefRHS(x1,dd); r-=2.0*staFunc->ScalProd(Pr,dd); staFunc->AddScaledVector(x,v0,-dddx/2,x1); h2_DefRHS(x1,dd); r+=2.0*staFunc->ScalProd(Pr,dd); staFunc->AddScaledVector(x,v0,-3*dddx/2,x1); h2_DefRHS(x1,dd); r-=(2.0/3.0)*staFunc->ScalProd(Pr,dd); staFunc->CopyVectorElements(Qi,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,3*dddx/2,x1); h2_DefRHS(x1,dd); r+=(2.0/3.0)*staFunc->ScalProd(Pi,dd); staFunc->AddScaledVector(x,v0,dddx/2,x1); h2_DefRHS(x1,dd); r-=2.0*staFunc->ScalProd(Pi,dd); staFunc->AddScaledVector(x,v0,-dddx/2,x1); h2_DefRHS(x1,dd); r+=2.0*staFunc->ScalProd(Pi,dd); staFunc->AddScaledVector(x,v0,-3*dddx/2,x1); h2_DefRHS(x1,dd); r-=(2.0/3.0)*staFunc->ScalProd(Pi,dd); staFunc->AddScaledVector(Qr,Qi,1.0,aa); staFunc->CopyVectorElements(aa,0,v0,0,nphase); staFunc->AddScaledVector(Pr,Pi,1.0,bb); staFunc->AddScaledVector(x,v0,3*dddx/2,x1); h2_DefRHS(x1,dd); r+=staFunc->ScalProd(bb,dd)/6; staFunc->AddScaledVector(x,v0,dddx/2,x1); h2_DefRHS(x1,dd); r-=staFunc->ScalProd(bb,dd)/2; staFunc->AddScaledVector(x,v0,-dddx/2,x1); h2_DefRHS(x1,dd); r+=staFunc->ScalProd(bb,dd)/2; staFunc->AddScaledVector(x,v0,-3*dddx/2,x1); h2_DefRHS(x1,dd); r-=staFunc->ScalProd(bb,dd)/6; staFunc->AddScaledVector(Qr,Qi,-1.0,aa); staFunc->CopyVectorElements(aa,0,v0,0,nphase); staFunc->AddScaledVector(Pr,Pi,-1.0,bb); staFunc->AddScaledVector(x,v0,3*dddx/2,x1); h2_DefRHS(x1,dd); r+=staFunc->ScalProd(bb,dd)/6; staFunc->AddScaledVector(x,v0,dddx/2,x1); h2_DefRHS(x1,dd); r-=staFunc->ScalProd(bb,dd)/2; staFunc->AddScaledVector(x,v0,-dddx/2,x1); h2_DefRHS(x1,dd); r+=staFunc->ScalProd(bb,dd)/2; staFunc->AddScaledVector(x,v0,-3*dddx/2,x1); h2_DefRHS(x1,dd); r-=staFunc->ScalProd(bb,dd)/6; r/=(dddx*dddx*dddx); } *res=(s+r)/(2*omega); } test[0]=*res; return 0; } Entry(Int2) h2_ProcBautin(Int2 Ind, Vector x, Vector v, CharPtr *msg) { stagenData->BifDataOutLen=0; if (Ind) sprintf(buf,"Generalized Hopf (?)"); else { /* reevaluate Jacobian */ if (JacRHS (x,C1)) return 1; staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->CopyVectorElements(VB1,0,V,0,nphase); staFunc->CopyVectorElements(VB2,0,W,0,nphase); /* orthogonalize and normalize eigenvectors */ staFunc->NormalizeVector(V); staFunc->AddScaledVector(W,V,-staFunc->ScalProd(V,W),L); staFunc->NormalizeVector(L); /* provide bifurcation data for gh -> h starter */ /* BifVector[0] to [nphase-1]: Hopf continuation vector*/ /* BifVector[nphase] to [2*nphase-1]: global vector L */ /* BifVector[2*nphase]: -lambda*lambda */ MemCopy(BifVector,V,sizeof(FloatHi)*nphase); MemCopy(BifVector+nphase,L,sizeof(FloatHi)*nphase); BifVector[2*nphase]=x[nphase+2+nappend]; stagenData->BifDataOutLen=sizeof(FloatHi)*(2*nphase+1); sprintf(buf,"Generalized Hopf: omega=%g",sqrt(BifVector[2*nphase])); } *msg=buf; return 0; } /* -lambda^2 */ Entry(Int2) h2_TestBogdanovTakens(Vector x, Vector v, FloatHiPtr res) { *res=x[nphase+2+nappend]; test[1]=*res; return 0; } Entry(Int2) h2_ProcBogdanovTakens(Int2 Ind, Vector x, Vector v, CharPtr *msg) { stagenData->BifDataOutLen=0; if (Ind) sprintf(buf,"Bogdanov-Takens (?)"); else { Vector v1,W0,W1; Matrix D; Int2 i,j,k; double a20,b20,b11,f0,f1,f2,f3,f4,f5,d; v1=staFunc->CreateVector(2*nphase); W0=staFunc->CreateVector(nphase); W1=staFunc->CreateVector(nphase); D=staFunc->CreateMatrix(2*nphase,2*nphase); /* compute (generalized)eigenvectors */ if (JacRHS (x,C1)) return 1; staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->CopyMatrixBlock(A,0,0,D,0,0,nphase,nphase); staFunc->CopyMatrixBlock(A,0,0,D,nphase,nphase,nphase,nphase); for (i=0; iSngSolve(D,_eps,v1); if (i != 2*nphase-2) { sprintf(buf,"Bogdanov-Takens: eigenvectors (?)"); goto final; } staFunc->CopyVectorElements(v1,0,V,0,nphase); /* null-vector */ d=sqrt(staFunc->ScalProd(V,V)); staFunc->NormalizeVector(V); staFunc->CopyVectorElements(v1,nphase,W,0,nphase); /* generalized null-vector */ staFunc->MultiplyVectorScalar(W,1.0/d,W); staFunc->AddScaledVector(W,V,-staFunc->ScalProd(V,W),W); /* provide bifurcation data for bt->* starters */ /* BifVector[0] to [nphase-1]: null-vector V (A*V=0, =1)*/ /* BifVector[nphase] to [2*nphase-1]: generalized null-vector W (A*W=V, =0) */ /* BifVector[2*nphase]: -lambda*lambda */ MemCopy(BifVector,V,sizeof(FloatHi)*nphase); MemCopy(BifVector+nphase,W,sizeof(FloatHi)*nphase); BifVector[2*nphase]=0.0; stagenData->BifDataOutLen=sizeof(FloatHi)*(2*nphase+1); /* compute the adjoint eigenvectors */ staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->TransposeMatrix(A,B); staFunc->ClearMatrix(D); staFunc->CopyMatrixBlock(B,0,0,D,0,0,nphase,nphase); staFunc->CopyMatrixBlock(B,0,0,D,nphase,nphase,nphase,nphase); for (i=0; iSngSolve(D,_eps,v1); if (i != 2*nphase-2) { sprintf(buf,"Bogdanov-Takens: adjoint eigenvectors (?)"); goto final; } staFunc->CopyVectorElements(v1,0,W1,0,nphase); /* adjoint null-vector W1 */ staFunc->CopyVectorElements(v1,nphase,W0,0,nphase); /* adjoint generalized null-vector W0 */ d=staFunc->ScalProd(V,W0); staFunc->MultiplyVectorScalar(W1,1.0/d,W1); staFunc->AddScaledVector(W0,W1,-staFunc->ScalProd(W0,W),W0); staFunc->MultiplyVectorScalar(W0,1.0/d,W0); /* compute the coefficients of the normal form : x'=y; y'=ax^2+bxy */ /* a20=, b20=, b11= */ if (stagenData->Der2) { CopyToX(x); CopyToP(x); H=stagenData->Der2(X,P,&zero); /* The 1st elem must be time! */ for (a20=0,b20=0,b11=0,i=0; iClearVector(v0); /* a20=, b20= */ staFunc->CopyVectorElements(V,0,v0,0,nphase); h2_DefRHS(x,F); f0=staFunc->ScalProd(W0,F); f3=staFunc->ScalProd(W1,F); staFunc->AddScaledVector(x,v0,ddx,x1); h2_DefRHS(x1,F); f1=staFunc->ScalProd(W0,F); f4=staFunc->ScalProd(W1,F); staFunc->AddScaledVector(x,v0,-ddx,x1); h2_DefRHS(x1,F); f2=staFunc->ScalProd(W0,F); f5=staFunc->ScalProd(W1,F); a20=(f1-2*f0+f2)/ddx/ddx; b20=(f4-2*f3+f5)/ddx/ddx; /* b11= */ staFunc->AddScaledVector(V,W,1.0,F); staFunc->CopyVectorElements(F,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,F); b11=staFunc->ScalProd(W1,F); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,F); b11+=staFunc->ScalProd(W1,F); staFunc->AddScaledVector(V,W,-1.0,F); staFunc->CopyVectorElements(F,0,v0,0,nphase); staFunc->AddScaledVector(x,v0,0.5*ddx,x1); h2_DefRHS(x1,F); b11-=staFunc->ScalProd(W1,F); staFunc->AddScaledVector(x,v0,-0.5*ddx,x1); h2_DefRHS(x1,F); b11-=staFunc->ScalProd(W1,F); b11/=(ddx*ddx); } if (b20>=0) sprintf(buf,"Bogdanov-Takens: a=%g, b=%g",b20/2,a20+b11); else sprintf(buf,"Bogdanov-Takens: a=%g, b=%g",-b20/2,-(a20+b11)); final: staFunc->DestroyMatrix(D); staFunc->DestroyVector(v1); W0=staFunc->DestroyVector(W0); W1=staFunc->DestroyVector(W1); } *msg=buf; return 0; } /* det of the Jacobian matrix */ Entry(Int2) h2_TestZeroHopf(Vector x, Vector v, FloatHiPtr res) { *res=1.0; if (nphase < 2) goto end; if (JacRHS (x,C1)) return 1; staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); if (staFunc->Decomp(A)) { *res=0.0; } else { staFunc->GetMatrixData(A,res,NULL); } end: test[2]=*res; return 0; } Entry(Int2) h2_ProcZeroHopf(Int2 Ind, Vector x, Vector v, CharPtr *msg) { stagenData->BifDataOutLen=0; if (Ind) sprintf(buf,"Zero-Hopf (?)"); else { /* reevaluate Jacobian */ if (JacRHS (x,C1)) return 1; staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->CopyVectorElements(VB1,0,V,0,nphase); staFunc->CopyVectorElements(VB2,0,W,0,nphase); /* orthogonalize and normalize eigenvectors */ staFunc->NormalizeVector(V); staFunc->AddScaledVector(W,V,-staFunc->ScalProd(V,W),L); staFunc->NormalizeVector(L); /* provide bifurcation data for zh -> h starter */ if (JacRHS (x,C1)) return 1; staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); if (staFunc->SngSolve(A,_eps,W) != nphase-1) { sprintf(buf,"Zero-Hopf: null-vector (?)"); *msg=buf; return 0; } staFunc->NormalizeVector(W); /* BifVector[0] to [nphase-1]: null-vector */ /* BifVector[nphase] to [2*nphase-1]: Hopf continuation vector */ /* BifVector[2*nphase] to [3*nphase-1]: global vector L */ /* BifVector[3*nphase]: -lambda*lambda */ MemCopy(BifVector,W,sizeof(FloatHi)*nphase); MemCopy(BifVector+nphase,V,sizeof(FloatHi)*nphase); MemCopy(BifVector+2*nphase,L,sizeof(FloatHi)*nphase); BifVector[3*nphase]=x[nphase+2+nappend]; stagenData->BifDataOutLen=sizeof(FloatHi)*(3*nphase+1); if (BifVector[3*nphase]>0) { sprintf(buf,"Zero-Hopf: omega=%g",sqrt(BifVector[3*nphase])); } else { sprintf(buf,"Zero eigenvalue + Neutral saddle: lambda=%g",sqrt(fabs(BifVector[3*nphase]))); } } *msg=buf; return 0; } /* DH test function */ Entry(Int2) h2_TestDoubleHopf(Vector x, Vector v, FloatHiPtr res) { VAR beta,lambda,phi,omega,d,s,r,maxV,maxW; Int2 i,ii,j,jj,imaxV,imaxW; *res=1.0; if (nphase < 4) goto end; if (JacRHS(x,C1)) { /* C1,A,AA,B,BB,V,W - global */ test[3]=*res; return 1; } if (x[nphase+2+nappend] <= 0.0) { /* neutral saddle */ lambda=sqrt(fabs(x[nphase+2+nappend])); /* compute real eigenvectors associated with lambda and -lambda */ staFunc->CopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); for (i=0; iSngSolve(AA,_eps,V) != nphase-1) { *res=1.0; test[3]=*res; return 1; } staFunc->CopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); for (i=0; iSngSolve(AA,_eps,W) != nphase-1) { *res=1.0; test[3]=*res; return 1; } staFunc->NormalizeVector(V); staFunc->AddScaledVector(W,V,-staFunc->ScalProd(V,W),W); /* ortogonalize */ staFunc->NormalizeVector(W); } else { /* Hopf */ omega=sqrt(x[nphase+2+nappend]); lambda=0.0; /* compute real and imaginary part of the adjoint eigenvector */ staFunc->CopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); staFunc->TransposeMatrix(AA,A); staFunc->ClearMatrix(BB); staFunc->CopyMatrixBlock(A,0,0,BB,0,0,nphase,nphase); staFunc->CopyMatrixBlock(A,0,0,BB,nphase,nphase,nphase,nphase); for (i=0; iSngSolve(BB,_eps,VV) != 2*nphase-2) { *res=1.0; test[3]=*res; return 1; } staFunc->CopyVectorElements(VV,0,V,0,nphase); staFunc->CopyVectorElements(VV,nphase,W,0,nphase); /* normalize real and imaginary parts =1, =1, =0 */ d=staFunc->ScalProd(V,V); s=staFunc->ScalProd(W,W); r=staFunc->ScalProd(V,W); if (r!=0) { beta=sqrt(fabs(4.0*r*r+(s-d)*(s-d))); if (2.0*fabs(2.0*r) > beta) { phi=PID2-atan(fabs((s-d)/(2.0*r))); } else { phi=atan(fabs(2.0*r/(s-d))); } if (r < 0) phi=PI-phi; if ((s-d) < 0) phi=-phi; phi/=2.0; for (i=0; iCopyVectorElements(VV,0,V,0,nphase); staFunc->CopyVectorElements(VV,nphase,W,0,nphase); } staFunc->NormalizeVector(V); staFunc->NormalizeVector(W); } /* find two basis vectors with max projection to Span(V,W) */ imaxV=0; maxV=fabs(V[0]); imaxW=0; maxW=fabs(W[0]); for (i=0; imaxV) { imaxV=i; maxV=fabs(V[i]); } if (fabs(W[i])>maxW) { imaxW=i; maxW=fabs(W[i]); } } if (imaxV==imaxW) { imaxW++; if (imaxW==nphase) imaxW=0; } /* fill in the matrix B with vectors to be orthogonalized as rows */ staFunc->ClearMatrix(B); staFunc->CopyVectorMatrixRow(V,B,0); staFunc->CopyVectorMatrixRow(W,B,1); for (i=0,ii=2,jj=0; iCopyMatrixRowVector(B,i,V); staFunc->CopyMatrixRowVector(B,i,F); for (j=0; jCopyMatrixRowVector(B,j,W); s=staFunc->ScalProd(V,W); staFunc->AddScaledVector(F,W,-s,F); } staFunc->NormalizeVector(F); staFunc->CopyVectorMatrixRow(F,B,i); } /* compute elements of the restriction matrix AH */ staFunc->CopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); staFunc->ClearMatrix(AH); for (i=2; iCopyMatrixRowVector(B,i,V); staFunc->MultiplyMatrixVector(AA,V,W); for (j=2; jCopyMatrixRowVector(B,j,V); AH[j-2][i-2]=staFunc->ScalProd(V,W); } } /* compute det of bi-product matrix DE */ staFunc->ClearMatrix(DE); { Int2 p,q,r,s; i=-1; for (p=1; pDecomp(DE)) { *res=0.0; } else { staFunc->GetMatrixData(DE,res,NULL); } end: test[3]=*res; return 0; } Entry(Int2) h2_ProcDoubleHopf(Int2 Ind, Vector x, Vector v, CharPtr *msg) { stagenData->BifDataOutLen=0; if (Ind) sprintf(buf,"Double Hopf (?)"); else {int i,j,imin,jmin; double r,s,d,phi,beta,gamma,lambda,omega; Vector Re,Im; Re=staFunc->CreateVector(nphase-2); Im=staFunc->CreateVector(nphase-2); staFunc->EigenValues(AH,Re,Im); /* find the pair of eigenvalues with zero sum */ s=fabs(Re[0]+Re[1])+fabs(Im[0]+Im[1]); imin=0; jmin=1; for (i=0; iCopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); for (i=0; iSngSolve(AA,_eps,V) != nphase-1) { sprintf(buf,"Double Hopf: lambda_2=%g, eigenvectors(?)",lambda); goto polufinal; } staFunc->CopyMatrixBlock(C1,0,0,AA,0,0,nphase,nphase); for (i=0; iSngSolve(AA,_eps,W) != nphase-1) { sprintf(buf,"Double Hopf: lambda_2=%g, eigenvectors(?)",lambda); goto polufinal; } staFunc->NormalizeVector(V); staFunc->NormalizeVector(W); staFunc->AddScaledVector(V,W,-1.0,F); staFunc->NormalizeVector(F); staFunc->AddScaledVector(V,W,1.0,V); staFunc->NormalizeVector(V); /* provide bifurcation data for dh->h starter */ /* Secondary branch: */ /* BifVector[2*nphase+1] to [3*nphase]: 2-Hopf continuation vector*/ /* BifVector[3*nphase+1] to [4*nphase]: 2-global vector L */ /* BifVector[4*nphase+1]: -lambda[2]*lambda[2] */ MemCopy(BifVector+2*nphase+1,V,sizeof(FloatHi)*nphase); MemCopy(BifVector+3*nphase+1,F,sizeof(FloatHi)*nphase); BifVector[4*nphase+1]=-lambda*lambda; stagenData->BifDataOutLen=sizeof(FloatHi)*(4*nphase+2); if (x[nphase+2+nappend] <= 0.0) { sprintf(buf, "Double neutral saddle: lambda[1]=%g, lambda[2]=%g",sqrt(fabs(x[nphase+2+nappend])),lambda); } else { sprintf(buf, "Hopf + neutral saddle: omega[1]=%g ,lambda_2=%g",sqrt(x[nphase+2+nappend]),lambda); } } else { /* Hopf */ lambda=Re[imin]; omega=fabs(Im[imin]); /* compute real and imaginary part of the eigenvector */ staFunc->ClearMatrix(BB); staFunc->CopyMatrixBlock(C1,0,0,BB,0,0,nphase,nphase); staFunc->CopyMatrixBlock(C1,0,0,BB,nphase,nphase,nphase,nphase); for (i=0; iSngSolve(BB,_eps,VV); if (i!= 2*nphase-2) { sprintf(buf,"Double Hopf: omega_2=%g, lambda_2=%g, eigenvectors(?)", omega,lambda); goto polufinal; } staFunc->CopyVectorElements(VV,0,V,0,nphase); staFunc->CopyVectorElements(VV,nphase,W,0,nphase); /* normalize real and imaginary parts */ d=staFunc->ScalProd(V,V); s=staFunc->ScalProd(W,W); r=staFunc->ScalProd(V,W); if (r!=0) { beta=sqrt(fabs(4.0*r*r+(s-d)*(s-d))); if (2.0*fabs(2.0*r) > beta) { phi=PID2-atan(fabs((s-d)/(2.0*r))); } else { phi=atan(fabs(2.0*r/(s-d))); } if (r < 0) phi=PI-phi; if ((s-d) < 0) phi=-phi; phi/=2.0; for (i=0; iCopyVectorElements(VV,0,V,0,nphase); d=staFunc->ScalProd(V,V); gamma=1/sqrt(d); staFunc->MultiplyVectorScalar(VV,gamma,VV); staFunc->CopyVectorElements(VV,nphase,F,0,nphase); staFunc->NormalizeVector(F); /* provide bifurcation data for dh->h starter */ /* Secondary branch: */ /* BifVector[2*nphase+1] to [3*nphase]: 2-Hopf continuation vector*/ /* BifVector[3*nphase+1] to [4*nphase]: 2-global vector L */ /* BifVector[4*nphase+1]: -lambda[2]*lambda[2] */ MemCopy(BifVector+2*nphase+1,VV,sizeof(FloatHi)*nphase); MemCopy(BifVector+3*nphase+1,F,sizeof(FloatHi)*nphase); BifVector[4*nphase+1]=omega*omega; stagenData->BifDataOutLen=sizeof(FloatHi)*(4*nphase+2); if (x[nphase+2+nappend] <= 0.0) { sprintf(buf, "Neutral saddle + Hopf: lambda[1]=%g, omega[2]=%g",sqrt(fabs(x[nphase+2+nappend])),omega); } else { sprintf(buf, "Double Hopf: omega[1]=%g, omega[2]=%g",sqrt(x[nphase+2+nappend]),omega); } } /* provide bifurcation data for dh -> h starter */ /* Primary branch: */ /* BifVector[0] to [nphase-1]: 1-Hopf continuation vector */ /* BifVector[nphase] to [2*nphase-1]: 1-global vector L */ /* BifVector[2*nphase]: -lambda[1]*lambda[1] */ MemCopy(BifVector,x+nphase+2+nappend,sizeof(FloatHi)*nphase); MemCopy(BifVector+nphase,L,sizeof(FloatHi)*nphase); BifVector[2*nphase]=x[nphase+2+nappend]; polufinal: Re=staFunc->DestroyVector(Re); Im=staFunc->DestroyVector(Im); } *msg=buf; return 0; } /**************************************************/ /* Adaptation of the bordering vectors V1,V2,W1,W2 */ Entry(Int2) h2_Adapter(Vector x, Vector v) { int i,j,l; VAR scal; Matrix S; /* update borders */ staFunc->TransposeMatrix(BP,BT); staFunc->Decomp(BT); staFunc->ClearVector(PP1); PP1[nphase]=1.0; staFunc->Solve(BT,PP1); staFunc->CopyVector(PP1,W1); staFunc->ClearVector(PP2); PP2[nphase+1]=1.0; staFunc->Solve(BT,PP2); staFunc->CopyVector(PP2,W2); staFunc->CopyVector(Q1,V1); staFunc->CopyVector(Q2,V2); staFunc->NormalizeVector(V1); staFunc->NormalizeVector(W1); staFunc->NormalizeVector(V2); staFunc->NormalizeVector(W2); /* update indices */ if (stagenData->Der2) { for (i=0;iDer2(X,P,&zero); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,VB,0,nphase); staFunc->MultiplyMatrixVector(A,VB,HV1); staFunc->CopyVectorElements(VB2,0,VB,0,nphase); staFunc->MultiplyMatrixVector(A,VB,HV2); staFunc->TransposeMatrix(A,B); staFunc->CopyVectorElements(PP1,0,VB,0,nphase); staFunc->MultiplyMatrixVector(B,VB,HW1); staFunc->CopyVectorElements(PP2,0,VB,0,nphase); staFunc->MultiplyMatrixVector(B,VB,HW2); S=staFunc->CreateMatrix(2,2); for (l=0; lDestroyMatrix(S); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(V1,0,HV1,0,nphase); staFunc->CopyVectorElements(V2,0,HV2,0,nphase); staFunc->CopyVectorElements(W1,0,HW1,0,nphase); staFunc->CopyVectorElements(W2,0,HW2,0,nphase); G0[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV1); G1[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV2); G2[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV1); G3[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV2); } else { staFunc->CopyMatrixBlock(C1,0,0,A,0,0,nphase,nphase); staFunc->ClearVector(HV1); staFunc->ClearVector(HV2); staFunc->ClearVector(HW1); staFunc->ClearVector(HW2); staFunc->CopyVectorElements(VB1,0,HV1,0,nphase); staFunc->CopyVectorElements(VB2,0,HV2,0,nphase); staFunc->MultiplyMatrixVector(A,HV1,HV1); staFunc->MultiplyMatrixVector(A,HV2,HV2); staFunc->TransposeMatrix(A,B); staFunc->CopyVectorElements(PP1,0,HW1,0,nphase); staFunc->CopyVectorElements(PP2,0,HW2,0,nphase); staFunc->MultiplyMatrixVector(B,HW1,HW1); staFunc->MultiplyMatrixVector(B,HW2,HW2); Jac1=staFunc->CreateMatrix(nphase,nphase); Jac2=staFunc->CreateMatrix(nphase,nphase); for (l=0; lClearVector(xp); staFunc->CopyVectorElements(x0,0,xp,0,nphase+2); xp[l]+=_dx; staFunc->ClearVector(xp1); staFunc->CopyVectorElements(x0,0,xp1,0,nphase+2); xp1[l]-=_dx; JacRHS(xp,C1); staFunc->CopyMatrixBlock(C1,0,0,Jac1,0,0,nphase,nphase); JacRHS(xp1,C1); staFunc->CopyMatrixBlock(C1,0,0,Jac2,0,0,nphase,nphase); /* g11 */ staFunc->CopyVectorElements(VB1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G0[l]=-staFunc->ScalProd(HW1,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G0[l]-=staFunc->ScalProd(K3,HV1); /* g12 */ staFunc->CopyVectorElements(VB2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G1[l]=-staFunc->ScalProd(HW1,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G1[l]-=staFunc->ScalProd(K3,HV2); /* g21 */ staFunc->CopyVectorElements(VB1,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G2[l]=-staFunc->ScalProd(HW2,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G2[l]-=staFunc->ScalProd(K3,HV1); /* g22 */ staFunc->CopyVectorElements(VB2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(Jac1,K3,K1); staFunc->MultiplyMatrixVector(Jac2,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G3[l]=-staFunc->ScalProd(HW2,K3); staFunc->TransposeMatrix(Jac1,B); staFunc->CopyVectorElements(PP2,0,K3,0,nphase); staFunc->MultiplyMatrixVector(B,K3,K1); staFunc->TransposeMatrix(Jac2,B); staFunc->MultiplyMatrixVector(B,K3,K2); staFunc->AddScaledVector(K1,K2,-1.0,K3); staFunc->MultiplyVectorScalar(K3,1.0/(_dx+_dx),K3); G3[l]-=staFunc->ScalProd(K3,HV2); } Jac1=staFunc->DestroyMatrix(Jac1); Jac2=staFunc->DestroyMatrix(Jac2); G0[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV1); G1[nphase+2+nappend]=-staFunc->ScalProd(HW1,HV2); G2[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV1); G3[nphase+2+nappend]=-staFunc->ScalProd(HW2,HV2); } for (i=0; iClearMatrix(JM2); JacRHS(x,C1); staFunc->CopyMatrixBlock(C1,0,0,JM2,0,0,nphase,nphase+2); for (i=0; iClearVector(TV2); TV2[nphase]=1.0; staFunc->Decomp(JM2); staFunc->Solve(JM2,TV2); staFunc->ClearMatrix(TVM); for (i=0;iClearVector(TV2); TV2[nphase+1]=1.0; staFunc->Solve(JM2,TV2); for (i=0;iClearVector(TV2); TV2[nphase+2]=1.0; staFunc->Solve(JM2,TV2); for (i=0;iTransposeMatrix(TVM,TVMT); staFunc->MultiplyMatrixMatrix(TVMT,TVM,HM); staFunc->Decomp(HM); staFunc->MultiplyMatrixVector(TVMT,G0,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G0); staFunc->MultiplyMatrixVector(TVMT,G1,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G1); staFunc->MultiplyMatrixVector(TVMT,G2,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G2); staFunc->MultiplyMatrixVector(TVMT,G3,HV); staFunc->Solve(HM,HV); staFunc->MultiplyMatrixVector(TVM,HV,G3); /* Determine component with maximum norm => gives m1,n1 */ G[0][0] = staFunc->NormOfVector(G0); G[0][1] = staFunc->NormOfVector(G1); G[1][0] = staFunc->NormOfVector(G2); G[1][1] = staFunc->NormOfVector(G3); m1=0; n1=0; if (G[0][1] > G[m1][n1]) n1=1; if (G[1][0] > G[m1][n1]) { m1=1; n1=0; } if (G[1][1] > G[m1][n1]) { m1=1; n1=1; } for (i=0;iNormalizeVector(G0); /* G1,G2,G3 are projected on G0 */ for (j=0;j<4;j++) { if (j != (2*m1+n1)) { for (i=0;iScalProd(G1,G0); staFunc->AddScaledVector(G1,G0,-scal,G1); for (i=0;i gives m2,n2 */ for (i=0;iNormOfVector(G0); G[0][1] = staFunc->NormOfVector(G1); G[1][0] = staFunc->NormOfVector(G2); G[1][1] = staFunc->NormOfVector(G3); G[m1][n1]=-1.0; m2=0; n2=0; if (G[0][1] > G[m2][n2]) n2=1; if (G[1][0] > G[m2][n2]) { m2=1; n2=0; } if (G[1][1] > G[m2][n2]) { m2=1; n2=1; } return 0; }