void computeHigherOrderChromaticities(LINE_LIST *beamline, double *clorb, RUN *run, long concatOrder, double deltaStep, long deltaPoints, long quickMode) { #define MAX_NDELTA_VALUES 101 /* must be at least 5 */ double trace[2][MAX_NDELTA_VALUES], delta[MAX_NDELTA_VALUES]; double eta[6]; double coef[MAX_NDELTA_VALUES], sCoef[MAX_NDELTA_VALUES], chi; long i, p; double c1; VMATRIX M1, M0, *Mp; if (deltaPoints>MAX_NDELTA_VALUES) bombElegant("too many points for higher-order chromaticity", NULL); if (deltaPoints<5) deltaPoints = 5; if (!(beamline->matrix)) bombElegant("no matrix for beamline (computeHigherOrderChromaticities)", NULL); beamline->chrom2[0] = beamline->chrom2[1] = 0; beamline->chrom3[0] = beamline->chrom3[1] = 0; initialize_matrices(&M0, 1); initialize_matrices(&M1, concatOrder); eta[0] = beamline->twiss0->etax; eta[1] = beamline->twiss0->etapx; eta[2] = beamline->twiss0->etay; eta[3] = beamline->twiss0->etapy; eta[4] = 0; eta[5] = 1; for (p=0; p<deltaPoints; p++) { delta[p] = (p-(deltaPoints/2+1))*deltaStep; for (i=0; i<6; i++) { M0.C[i] = (clorb?clorb[i]:0)+delta[p]*eta[i] + (i<4? sqr(delta[p])*beamline->eta2[i] + pow3(delta[p])*beamline->eta3[i] : 0); M0.R[i][i] = 1; } if (quickMode) concat_matrices(Mp=&M1, beamline->matrix, &M0, 0); else Mp = append_full_matrix(beamline->elem_twiss, run, &M0, concatOrder); for (i=0; i<2; i++) /* Tr[0,1][p] is the trace for x,y plane for point p */ trace[i][p] = Mp->R[2*i][2*i] + Mp->R[2*i+1][2*i+1]; if (!quickMode) { free_matrices(Mp); free(Mp); Mp = NULL; } } for (i=0; i<2; i++) { lsfn(delta, trace[i], NULL, deltaPoints, (long)(MIN(deltaPoints-2,5)), coef, sCoef, &chi, NULL); c1 = -coef[1]/(4*PI*sin(PIx2*beamline->tune[i])); beamline->chrom2[i] = ( 2*coef[2] + 8*sqr(PI)*sqr(c1)*cos(PIx2*beamline->tune[i]) )/(-4*PI*sin(PIx2*beamline->tune[i])); if (quickMode) beamline->chrom3[i] = 0; else beamline->chrom3[i] = ( 6*coef[3] - 16*pow3(PI*c1)*sin(PIx2*beamline->tune[i]) + 24*sqr(PI)*c1*beamline->chrom2[i]*cos(PIx2*beamline->tune[i]) )/(-4*PI*sin(PIx2*beamline->tune[i])); } free_matrices(&M0); free_matrices(&M1); }
/* free all of our data structures */ void cleanup(system_t *system) { int i,j; #ifdef QM_ROTATION if(system->quantum_rotation) free_rotational(system); #endif /* QM_ROTATION */ if(system->polarization && !system->cuda) free_matrices(system); if(system->polar_wolf_alpha_lookup) if ( system->polar_wolf_alpha_table ) free(system->polar_wolf_alpha_table); //need to rebuild atom and pair arrays so we can free everything system->natoms = countNatoms(system); rebuild_arrays(system); free_all_pairs(system); free_all_molecules(system, system->molecules); //free our arrays free(system->molecule_array); free(system->atom_array); free(system->pqr_input); free(system->pqr_output); free(system->energy_output); free(system->energy_output_csv); free(system->virial_output); if(system->surf_output) free(system->surf_output); if(system->traj_output) free(system->traj_output); if(system->traj_input) free(system->traj_input); if(system->dipole_output) free(system->dipole_output); if(system->field_output) free(system->field_output); if(system->frozen_output) free(system->frozen_output); if(system->surf_preserve_rotation_on) free(system->surf_preserve_rotation_on); if(system->cavity_bias) free_cavity_grid(system); if ( system->vdw_eiso_info ) free_vdw_eiso(system->vdw_eiso_info); // free multi sorbate related stuff if ( system->fugacities ) free(system->fugacities); if(system->insert_input) free(system->insert_input); //insert.pdb arrays and shit if(system->insertion_molecules_array) free(system->insertion_molecules_array); if(system->insertion_molecules) free_all_molecules(system, system->insertion_molecules); // free sorbate info array free(system->sorbateInfo); /* free statistics */ free(system->nodestats); free(system->avg_nodestats); free(system->observables); free(system->avg_observables); free(system->checkpoint->observables); if ( system->checkpoint->molecule_backup != NULL ) free_molecule(system, system->checkpoint->molecule_backup); free(system->checkpoint); /*free histogram stuff*/ for ( i=0; i<system->grids->histogram->x_dim; i++ ) { for ( j=0; j<system->grids->histogram->y_dim; j++ ) { free(system->grids->histogram->grid[i][j]); } free(system->grids->histogram->grid[i]); } if ( !rank ) { for ( i=0; i<system->grids->histogram->x_dim; i++ ) { for ( j=0; j<system->grids->histogram->y_dim; j++ ) { free(system->grids->avg_histogram->grid[i][j]); } free(system->grids->avg_histogram->grid[i]); } free(system->grids->avg_histogram->grid); } free(system->grids->histogram->grid); free(system->grids->avg_histogram); free(system->grids->histogram); free(system->grids); free(system->histogram_output); /* free fit input lists*/ if ( system->fit_input_list.data.count > 0 ) { fileNode_t *node = &(system->fit_input_list); fileNode_t *nextnode; //the first one is statically declared //the others are malloc'd in input.c if ( node->next ) node=node->next; while ( node->next ) { nextnode=node->next; free(node->data.filename); free(node); node=nextnode; } free(node->data.filename); free(node); } free(system->pqr_restart); free(system->pbc); free(system->job_name); // (CRC) kill_rng(); free(system); }
long do_chromaticity_correction(CHROM_CORRECTION *chrom, RUN *run, LINE_LIST *beamline, double *clorb, long step, long last_iteration) { VMATRIX *M; double chromx0, chromy0, dchromx=0, dchromy=0; double K2=0.0, *K2ptr; ELEMENT_LIST *context; long i, K2_param=0, type=0, iter, count; double beta_x, alpha_x, eta_x, etap_x; double beta_y, alpha_y, eta_y, etap_y; double K2_min, K2_max; unsigned long unstable; double lastError, presentError; char buffer[256]; log_entry("do_chromaticity_correction"); if (!beamline->matrix || !beamline->twiss0) { if (!beamline->twiss0) beamline->twiss0 = tmalloc(sizeof(*beamline->twiss0)); if (!beamline->elem_twiss) { ELEMENT_LIST *eptr; eptr = beamline->elem_twiss = &(beamline->elem); while (eptr) { if (eptr->type==T_RECIRC) beamline->elem_twiss = beamline->elem_recirc = eptr; eptr = eptr->succ; } } if (beamline->matrix) { free_matrices(beamline->matrix); free(beamline->matrix); beamline->matrix = NULL; } beamline->matrix = compute_periodic_twiss(&beta_x, &alpha_x, &eta_x, &etap_x, beamline->tune, &beta_y, &alpha_y, &eta_y, &etap_y, beamline->tune+1, beamline->elem_twiss, clorb, run, &unstable, NULL, NULL); beamline->twiss0->betax = beta_x; beamline->twiss0->alphax = alpha_x; beamline->twiss0->phix = 0; beamline->twiss0->etax = eta_x; beamline->twiss0->etapx = etap_x; beamline->twiss0->betay = beta_y; beamline->twiss0->alphay = alpha_y; beamline->twiss0->phiy = 0; beamline->twiss0->etay = eta_y; beamline->twiss0->etapy = etap_y; propagate_twiss_parameters(beamline->twiss0, beamline->tune, beamline->waists, NULL, beamline->elem_twiss, run, clorb, beamline->couplingFactor); } else if (beamline->matrix->order<2) { if (beamline->matrix) { free_matrices(beamline->matrix); free(beamline->matrix); beamline->matrix = NULL; } beamline->matrix = full_matrix(beamline->elem_twiss, run, 2); } if (!(M = beamline->matrix) || !M->C || !M->R || !M->T) bombElegant("something wrong with transfer map for beamline (do_chromaticity_correction.1)", NULL); computeChromaticities(&chromx0, &chromy0, NULL, NULL, NULL, NULL, beamline->twiss0, beamline->elast->twiss, M); if (verbosityLevel>0) { fprintf(stdout, "\nAdjusting chromaticities:\n"); fflush(stdout); fprintf(stdout, "initial chromaticities: %e %e\n", chromx0, chromy0); fflush(stdout); } presentError = DBL_MAX; for (iter=0; iter<chrom->n_iterations; iter++) { K2_max = -(K2_min = DBL_MAX); dchromx = chrom->chromx - chromx0; dchromy = chrom->chromy - chromy0; if (chrom->tolerance>0 && chrom->tolerance>fabs(dchromx) && chrom->tolerance>fabs(dchromy) ) break; lastError = presentError; presentError = sqr(dchromx)+sqr(dchromy); if (iter && presentError>lastError) { fprintf(stdout, "Error increasing---reducing gain\n"); fflush(stdout); chrom->correction_fraction /= 10; } if (chrom->use_perturbed_matrix) computeChromCorrectionMatrix(run, beamline, chrom); chrom->dchrom->a[0][0] = dchromx; chrom->dchrom->a[1][0] = dchromy; m_mult(chrom->dK2, chrom->T, chrom->dchrom); for (i=0; i<chrom->n_families; i++) { if (isnan(chrom->correction_fraction*chrom->dK2->a[i][0]) || isinf(chrom->correction_fraction*chrom->dK2->a[i][0])) break; } if (i!=chrom->n_families) { fprintf(stdout, "Unable to correct chromaticity---diverged\n"); fflush(stdout); return 0; } for (i=0; i<chrom->n_families; i++) { context = NULL; count = 0; while ((context=find_element(chrom->name[i], &context, beamline->elem_twiss))) { if (count==0 && (K2_param = confirm_parameter("K2", context->type))<0) { fprintf(stdout, "error: element %s doesn't have K2 parameter\n", context->name); fflush(stdout); exitElegant(1); } if (!(K2ptr = (double*)(context->p_elem + entity_description[context->type].parameter[K2_param].offset))) bombElegant("K2ptr NULL in setup_chromaticity_correction", NULL); K2 = (*K2ptr += chrom->correction_fraction*chrom->dK2->a[i][0]); if (chrom->strengthLimit>0 && chrom->strengthLimit<fabs(K2)) { K2 = *K2ptr = SIGN(K2)*chrom->strengthLimit; } sprintf(buffer, "%s#%ld.K2", context->name, context->occurence); rpn_store(K2, NULL, rpn_create_mem(buffer, 0)); if (K2>K2_max) K2_max = K2; if (K2<K2_min) K2_min = K2; if (context->matrix) { free_matrices(context->matrix); free(context->matrix); context->matrix = NULL; } compute_matrix(context, run, NULL); type = context->type; count++; /* fprintf(stdout, "new value of %s#%ld[K2] is %.15lg 1/m^3\n", chrom->name[i], context->occurence, K2); fflush(stdout); */ } if (verbosityLevel>1) { fprintf(stdout, "Change for family %ld (%ld sextupoles): %e\n", i, count, chrom->correction_fraction*chrom->dK2->a[i][0]); fflush(stdout); } if (alter_defined_values) change_defined_parameter(chrom->name[i], K2_param, type, K2, NULL, LOAD_FLAG_ABSOLUTE); } if (beamline->links) { /* rebaseline_element_links(beamline->links, run, beamline); */ assert_element_links(beamline->links, run, beamline, STATIC_LINK+DYNAMIC_LINK+(alter_defined_values?LINK_ELEMENT_DEFINITION:0)); } if (beamline->matrix) { free_matrices(beamline->matrix); free(beamline->matrix); beamline->matrix = NULL; } M = beamline->matrix = compute_periodic_twiss(&beta_x, &alpha_x, &eta_x, &etap_x, beamline->tune, &beta_y, &alpha_y, &eta_y, &etap_y, beamline->tune+1, beamline->elem_twiss, clorb, run, &unstable, NULL, NULL); beamline->twiss0->betax = beta_x; beamline->twiss0->alphax = alpha_x; beamline->twiss0->phix = 0; beamline->twiss0->etax = eta_x; beamline->twiss0->etapx = etap_x; beamline->twiss0->betay = beta_y; beamline->twiss0->alphay = alpha_y; beamline->twiss0->phiy = 0; beamline->twiss0->etay = eta_y; beamline->twiss0->etapy = etap_y; if (!M || !M->C || !M->R || !M->T) bombElegant("something wrong with transfer map for beamline (do_chromaticity_correction.2)", NULL); computeChromaticities(&chromx0, &chromy0, NULL, NULL, NULL, NULL, beamline->twiss0, beamline->elast->twiss, M); beamline->chromaticity[0] = chromx0; beamline->chromaticity[1] = chromy0; if (verbosityLevel>0) { fprintf(stdout, "resulting chromaticities: %e %e\n", chromx0, chromy0); fprintf(stdout, "min, max sextupole strength: %e %e 1/m^2\n", K2_min, K2_max); fflush(stdout); } } if (fp_sl && last_iteration) { for (i=0; i<chrom->n_families; i++) { context = NULL; while ((context=find_element(chrom->name[i], &context, beamline->elem_twiss))) { if (( K2_param = confirm_parameter("K2", context->type))<0) bombElegant("confirm_parameter doesn't return offset for K2 parameter.\n", NULL); fprintf(fp_sl, "%ld %e %s\n", step, *((double*)(context->p_elem + entity_description[context->type].parameter[K2_param].offset)), chrom->name[i]); } } fflush(fp_sl); } propagate_twiss_parameters(beamline->twiss0, beamline->tune, beamline->waists, NULL, beamline->elem_twiss, run, clorb, beamline->couplingFactor); log_exit("do_chromaticity_correction"); if (chrom->tolerance>0 && (chrom->tolerance<fabs(dchromx) || chrom->tolerance<fabs(dchromy)) && chrom->exit_on_failure) { fprintf(stdout, "Chromaticity correction failure---exiting!\n"); exitElegant(1); } return 1; }
void computeChromCorrectionMatrix(RUN *run, LINE_LIST *beamline, CHROM_CORRECTION *chrom) { VMATRIX *M; double chromx, chromy; double chromx0, chromy0; double K2=0.0, *K2ptr; ELEMENT_LIST *context; long i, count, K2_param=0; MATRIX *C, *Ct, *CtC, *inv_CtC; m_alloc(&C, 2, chrom->n_families); m_alloc(&Ct, chrom->n_families, 2); m_alloc(&CtC, chrom->n_families, chrom->n_families); m_alloc(&inv_CtC, chrom->n_families, chrom->n_families); if (chrom->T) m_free(&(chrom->T)); if (chrom->dK2) m_free(&(chrom->dK2)); if (chrom->dchrom) m_free(&(chrom->dchrom)); m_alloc(&(chrom->T), chrom->n_families, 2); m_alloc(&(chrom->dK2), chrom->n_families, 1); m_alloc(&(chrom->dchrom), 2, 1); if (verbosityLevel>2) { fprintf(stdout, "Computing chromaticity influence matrix for all named sextupoles.\n"); fflush(stdout); } computeChromaticities(&chromx0, &chromy0, NULL, NULL, NULL, NULL, beamline->twiss0, beamline->elast->twiss, M=beamline->matrix); M = NULL; for (i=0; i<chrom->n_families; i++) { count = 0; context = NULL; while ((context=find_element(chrom->name[i], &context, beamline->elem_twiss))) { if (!count && !(K2_param=confirm_parameter("K2", context->type))) { fprintf(stdout, "error: element %s does not have K2 parameter\n", context->name); fflush(stdout); exitElegant(1); } if (!(K2ptr = (double*)(context->p_elem + entity_description[context->type].parameter[K2_param].offset))) bombElegant("K2ptr NULL in setup_chromaticity_correction", NULL); if (count==0) K2 = *K2ptr; *K2ptr = K2 + chrom->sextupole_tweek; if (context->matrix) { free_matrices(context->matrix); free(context->matrix); context->matrix = NULL; } compute_matrix(context, run, NULL); count++; } if (count==0) { fprintf(stdout, "error: element %s is not in the beamline.\n", chrom->name[i]); fflush(stdout); exitElegant(1); } if (beamline->links) { /* rebaseline_element_links(beamline->links, run, beamline); */ assert_element_links(beamline->links, run, beamline, STATIC_LINK+DYNAMIC_LINK); } if (M) { free_matrices(M); free(M); M = NULL; } M = full_matrix(beamline->elem_twiss, run, 2); computeChromaticities(&chromx, &chromy, NULL, NULL, NULL, NULL, beamline->twiss0, beamline->elast->twiss, M); C->a[0][i] = (chromx-chromx0)/chrom->sextupole_tweek; C->a[1][i] = (chromy-chromy0)/chrom->sextupole_tweek; if (C->a[0][i]==0 || C->a[1][i]==0) { fprintf(stdout, "error: element %s does not change the chromaticity!\n", chrom->name[i]); fflush(stdout); exitElegant(1); } count = 0; context = NULL; while ((context=find_element(chrom->name[i], &context, beamline->elem_twiss))) { if (!count && !(K2_param=confirm_parameter("K2", context->type))) { fprintf(stdout, "error: element %s does not have K2 parameter\n", context->name); fflush(stdout); exitElegant(1); } if (!(K2ptr = (double*)(context->p_elem + entity_description[context->type].parameter[K2_param].offset))) bombElegant("K2ptr NULL in setup_chromaticity_correction", NULL); if (!K2ptr) bombElegant("K2ptr NULL in setup_chromaticity_correction", NULL); *K2ptr = K2; if (context->matrix) { free_matrices(context->matrix); free(context->matrix); context->matrix = NULL; } compute_matrix(context, run, NULL); count++; } } if (M) { free_matrices(M); free(M); M = NULL; } if (beamline->matrix) { free_matrices(beamline->matrix); free(beamline->matrix); beamline->matrix = NULL; } beamline->matrix = full_matrix(beamline->elem_twiss, run, run->default_order); if (verbosityLevel>1) { fprintf(stdout, "\nfamily dCHROMx/dK2 dCHROMy/dK2\n"); fflush(stdout); for (i=0; i<chrom->n_families; i++) fprintf(stdout, "%10s: %14.7e %14.7e\n", chrom->name[i], C->a[0][i], C->a[1][i]); fflush(stdout); } m_trans(Ct, C); m_mult(CtC, Ct, C); m_invert(inv_CtC, CtC); m_mult(chrom->T, inv_CtC, Ct); if (verbosityLevel>1) { fprintf(stdout, "\nfamily dK2/dCHROMx dK2/dCHROMy\n"); fflush(stdout); for (i=0; i<chrom->n_families; i++) fprintf(stdout, "%10s: %14.7e %14.7e\n", chrom->name[i], chrom->T->a[i][0], chrom->T->a[i][1]); fflush(stdout); fprintf(stdout, "\n"); fflush(stdout); } m_free(&C); m_free(&Ct); m_free(&CtC); m_free(&inv_CtC); }
int run_coupled_twiss_output(RUN *run, LINE_LIST *beamline, double *starting_coord) { char JOBVL, JOBVR; int N, LDA, LDVL, LDVR, lwork, info, i, j, k; double A[36], WR[6], WI[6], VL[36], VR[36], work[1000]; double emit[3], Norm[3], Vnorm[36]; double Amatrix[108], SigmaMatrix[6][6]; int matDim, eigenModesNumber; double transferMatrix[36]; VMATRIX *M, *M1; double **R; ELEMENT_LIST *eptr, *eptr0; long nElements, lastNElements, iElement; double betax1, betax2, betay1, betay2, etax, etay, tilt; if (!initialized) return 0; if (verbosity>1) fprintf(stdout, "\n* Computing coupled sigma matrix\n"); if (emittances_from_twiss_command) { if (!(beamline->flags&BEAMLINE_TWISS_DONE)) { fprintf(stderr, "emittances_from_twiss_command was set but twiss calculations not seen"); return(1); } if (!(beamline->flags&BEAMLINE_RADINT_DONE)) { fprintf(stderr, "emittances_from_twiss_command was set but radiation integral calculations not seen"); return(1); } emit_x = beamline->radIntegrals.ex0; sigma_dp = beamline->radIntegrals.sigmadelta; if (verbosity>1) fprintf(stdout, "Raw emittance = %e, momentum spread = %e\n", emit_x, sigma_dp); } fflush(stdout); emit[0] = emit_x; emit[1] = emit_x*emittance_ratio; /* Count the number of elements from the recirc element to the end. */ /* Also store the pointer to the recirc element. */ eptr = eptr0 = &(beamline->elem); nElements = lastNElements = beamline->n_elems; while (eptr) { if (eptr->type==T_RECIRC) { lastNElements = nElements; eptr0 = eptr; } eptr = eptr->succ; nElements--; } nElements = lastNElements; if (starting_coord) { /* use the closed orbit to compute the on-orbit R matrix */ M1 = tmalloc(sizeof(*M1)); initialize_matrices(M1, 1); for (i=0; i<6; i++) { M1->C[i] = starting_coord[i]; M1->R[i][i] = 1; } M = accumulate_matrices(eptr0, run, M1, concat_order, 0); free_matrices(M1); free(M1); M1 = NULL; } else M = accumulate_matrices(eptr0, run, NULL, concat_order, 0); R = M->R; if (verbosity > 2) { long order; order = M->order; M->order = 1; print_matrices(stdout, "One-turn matrix:", M); M->order = order; } /* Determination of matrix dimension for these calculations. */ if (calculate_3d_coupling != 1) { matDim=4; } else { if (abs(R[4][4])+abs(R[5][5])>=2) { printf("Either there is no cavity or 3rd mode is unstable. Only 2 modes will be calculated.\n"); matDim=4; } else { matDim=6; } } eigenModesNumber=matDim/2; /*--- Reducing matrix dimensions, A is reduced R */ for (i=0; i<matDim; i++) { for (j=0; j<matDim; j++) { A[i*matDim+j]=R[j][i]; } } free_matrices(M); free(M); M = NULL; /*--- Changing time sign for symplecticity... */ if (matDim == 6) { for (i=0; i<6; i++) { A[24+i]=-1.0*A[24+i]; A[i*6+4]=-1.0*A[i*6+4]; } } if (verbosity > 3) { MatrixPrintout((double*)&A, &matDim, &matDim, 1); } /*--- Calculating eigenvectors using dgeev_ ... */ JOBVL='N'; JOBVR='V'; N=matDim; LDA=matDim; LDVL=1; LDVR=matDim; lwork=204; #if defined(SUNPERF) || defined(LAPACK) || defined(CLAPACK) dgeev_((char*)&JOBVL, (char*)&JOBVR, (int*)&N, (double*)&A, (int*)&LDA, (double*)&WR, (double*)&WI, (double*)&VL, (int*)&LDVL, (double*)&VR, (int*)&LDVR, (double*)&work, (int*)&lwork, (int*)&info); #else fprintf(stderr, "Error calling dgeev. You will need to install LAPACK and rebuild elegant\n"); return(1); #endif if (info != 0) { if (info < 0) { printf("Error calling dgeev, argument %d.\n", abs(info)); } if (info > 0) { printf("Error running dgeev, calculation of eigenvalue number %d failed.\n", info); } return(1); } if (verbosity > 0) { printf("Info: %d ; %f \n", info, work[0]); for(i=0; i<matDim; i++) { printf("%d: %9.6f + i* %10.6f\n",i,WR[i],WI[i]); } fflush(stdout); } if (verbosity > 1) { printf("Non-normalized vectors:\n"); MatrixPrintout((double*)&VR, &matDim, &matDim, 1); fflush(stdout); } /*--- Sorting of eigenvalues and eigenvectors according to (x,y,z)... */ SortEigenvalues((double*)&WR, (double*)&WI, (double*)&VR, matDim, eigenModesNumber, verbosity); /*--- Normalization of eigenvectors... */ for (k=0; k<eigenModesNumber; k++) { Norm[k]=0; for (i=0; i<eigenModesNumber; i++) { /* Index = Irow*matDim + Icolumn */ Norm[k]+=VR[2*k*matDim+2*i+1]*VR[(2*k+1)*matDim+2*i]-VR[2*k*matDim+2*i]*VR[(2*k+1)*matDim+2*i+1]; } Norm[k]=1.0/sqrt(fabs(Norm[k])); if (verbosity > 2) { printf("Norm[%d]= %12.4e \n",k,Norm[k]); } } for (k=0; k<eigenModesNumber; k++) { for (i=0; i<matDim; i++) { Vnorm[k*2*matDim+i]=VR[k*2*matDim+i]*Norm[k]; Vnorm[(k*2+1)*matDim+i]=VR[(k*2+1)*matDim+i]*Norm[k]; } } if (verbosity > 1) { printf("Normalized vectors:\n"); MatrixPrintout((double*)&Vnorm, &matDim, &matDim, 1); } if (SDDScoupledInitialized) { /*--- Prepare the output file */ if (!SDDS_StartPage(&SDDScoupled, nElements)) { fflush(stdout); SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors); return(1); } } /*--- Loop over elements */ iElement=0; eptr = eptr0; while (eptr) { if (verbosity > 0) { printf("\nElement number %ld: %s\n", iElement, eptr->name); fflush(stdout); } if (!eptr->accumMatrix) { fprintf(stderr, "Error: no accumulated matrix found for element %s", eptr->name); return(1); } /*--- Reducing matrix dimensions */ R = eptr->accumMatrix->R; for (i=0; i<matDim; i++) { for (j=0; j<matDim; j++) { transferMatrix[i*matDim+j]=R[j][i]; } } /*--- Changing time sign for symplecticity... */ if (matDim == 6) { for (i=0; i<6; i++) { transferMatrix[24+i]= -1.0*transferMatrix[24+i]; transferMatrix[i*6+4]=-1.0*transferMatrix[i*6+4]; } } /*--- Calculating A matrices (product of eigenvectors)... */ GetAMatrix((double*)&Vnorm, (double*)&transferMatrix, (double*)&Amatrix, &eigenModesNumber, &matDim); if (verbosity > 1) { for (k=0; k<eigenModesNumber; k++) { printf("A matrix for mode %d\n", k); MatrixPrintout((double*)&Amatrix[k*matDim*matDim], &matDim, &matDim, 1); } } /*--- Calculating sigma matrix... */ if (eigenModesNumber == 3) { emit[2]=sigma_dp*sigma_dp*Amatrix[2*matDim*matDim+4*matDim+4]; } for (i=0; i<matDim; i++) { for (j=0; j<matDim; j++) { SigmaMatrix[i][j]=0; for (k=0; k<eigenModesNumber; k++) { SigmaMatrix[i][j]+=emit[k]*Amatrix[k*matDim*matDim+i*matDim+j]; } } } if (verbosity > 0) { printf("Sigma matrix:\n"); MatrixPrintout((double*)&SigmaMatrix, &matDim, &matDim, 2); } tilt=0.5*atan(2*SigmaMatrix[0][2]/(SigmaMatrix[0][0]-SigmaMatrix[2][2])); if (SDDScoupledInitialized) { /*--- Calculating beam sizes: 0-SigmaX, 1-SigmaXP, 2-SigmaY, 3-SigmaYP, 4-BeamTilt, 5-BunchLength */ if (!SDDS_SetRowValues(&SDDScoupled, SDDS_SET_BY_NAME|SDDS_PASS_BY_VALUE, iElement, "ElementName", eptr->name, "s", eptr->end_pos, "Sx", sqrt(SigmaMatrix[0][0]), "Sxp", sqrt(SigmaMatrix[1][1]), "Sy", sqrt(SigmaMatrix[2][2]), "Syp", sqrt(SigmaMatrix[3][3]), "xyTilt", tilt, "Ss", eigenModesNumber==3?sqrt(SigmaMatrix[4][4]):-1, NULL)) { SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors); return(1); } } if (verbosity > 0) { printf("SigmaX = %12.4e, SigmaY = %12.4e, Beam tilt = %12.4e \n", sqrt(SigmaMatrix[0][0]), sqrt(SigmaMatrix[2][2]), 0.5*atan(2*SigmaMatrix[0][2]/(SigmaMatrix[0][0]-SigmaMatrix[2][2]))); printf("SigmaXP = %12.4e, SigmaYP = %12.4e, \n", sqrt(SigmaMatrix[1][1]), sqrt(SigmaMatrix[3][3])); if (eigenModesNumber==3) { printf("Bunch length = %12.4e \n", sqrt(SigmaMatrix[4][4])); } } betax1 = Amatrix[0]; betax2 = Amatrix[1*matDim*matDim]; betay1 = Amatrix[2*matDim+2]; betay2 = Amatrix[1*matDim*matDim+2*matDim+2]; etax = sqrt(Amatrix[2*matDim*matDim]*Amatrix[2*matDim*matDim+4*matDim+4]); etay = sqrt(Amatrix[2*matDim*matDim+2*matDim+2]*Amatrix[2*matDim*matDim+4*matDim+4]); if (SDDScoupledInitialized) { if (!SDDS_SetRowValues(&SDDScoupled, SDDS_SET_BY_NAME|SDDS_PASS_BY_VALUE, iElement, "betax1", betax1, "betax2", betax2, "betay1", betay1, "betay2", betay2, "etax", etax, "etay", etay, NULL)) { SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors); return(1); } if (output_sigma_matrix) { char name[100]; for (i=0; i<matDim; i++) for (j=i; j<matDim; j++) { sprintf(name, "S%d%d", i+1, j+1); if (!SDDS_SetRowValues(&SDDScoupled, SDDS_SET_BY_NAME|SDDS_PASS_BY_VALUE, iElement, name, SigmaMatrix[i][j], NULL)) { SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors); return(1); } } } } if (verbosity > 0) { printf("betax_1 = %12.4e, betax_2 = %12.4e \n", Amatrix[0], Amatrix[1*matDim*matDim]); printf("betay_1 = %12.4e, betay_2 = %12.4e \n", Amatrix[2*matDim+2], Amatrix[1*matDim*matDim+2*matDim+2]); printf("etax = %12.4e, etay = %12.4e \n", sqrt(Amatrix[2*matDim*matDim]*Amatrix[2*matDim*matDim+4*matDim+4]), sqrt(Amatrix[2*matDim*matDim+2*matDim+2]*Amatrix[2*matDim*matDim+4*matDim+4])); fflush(stdout); } if (eptr->type==T_MARK && ((MARK*)eptr->p_elem)->fitpoint) store_fitpoint_ctwiss_parameters((MARK*)eptr->p_elem, eptr->name, eptr->occurence, betax1, betax2, betay1, betay2, etax, etay, tilt); iElement++; eptr = eptr->succ; } if (SDDScoupledInitialized && !SDDS_WritePage(&SDDScoupled)) { SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors); return(1); } return(0); }
long assert_element_links(ELEMENT_LINKS *links, RUN *run, LINE_LIST *beamline, long flags) { long i_link, i_elem, i_item, matrices_changed; long elem_type, data_type, param; double value; ELEMENT_LIST **targ, **sour; char *p_elem; short lengthChanged = 0; log_entry("assert_element_links"); if (!links || links->n_links==0) { log_exit("assert_element_links"); return(0); } if (!links->target_name || !links->item || !links->source_name || !links->equation || !links->n_targets || !links->target_elem || !links->source_elem) { fputs("error: link structure has null pointers (assert_element_links)", stdout); abort(); } matrices_changed = 0; for (i_link=0; i_link<links->n_links; i_link++) { if (!(flags&links->flags[i_link])) continue; targ = links->target_elem[i_link]; sour = links->source_elem[i_link]; elem_type = targ[0]->type; param = links->target_param[i_link]; data_type = entity_description[elem_type].parameter[param].type; #if DEBUG fprintf(stdout, "asserting %ld links of %s.%s to %s\n", links->n_targets[i_link], links->target_name[i_link], links->item[i_link], links->source_name[i_link]); fflush(stdout); fprintf(stdout, "source type is %ld, with %ld parameters\n", sour[0]->type, entity_description[sour[0]->type].n_params); fflush(stdout); #endif for (i_elem=0; i_elem<links->n_targets[i_link]; i_elem++) { #if DEBUG fprintf(stdout, " working on element %ld\n", i_elem); fflush(stdout); #endif p_elem = sour[i_elem]->p_elem; for (i_item=0; i_item<entity_description[sour[0]->type].n_params; i_item++) { char s[1024]; double value0; switch (entity_description[sour[0]->type].parameter[i_item].type) { case IS_DOUBLE: value = *((double*)(p_elem+entity_description[sour[0]->type].parameter[i_item].offset)); rpn_store(value, NULL, rpn_create_mem(entity_description[sour[0]->type].parameter[i_item].name, 0)); value0 = *((double*)(sour[i_elem]->p_elem0+entity_description[sour[0]->type].parameter[i_item].offset)); sprintf(s, "%s0", entity_description[sour[0]->type].parameter[i_item].name); rpn_store(value0, NULL, rpn_create_mem(s, 0)); break; case IS_LONG: value = *((long*)(p_elem+entity_description[sour[0]->type].parameter[i_item].offset)); rpn_store(value, NULL, rpn_create_mem(entity_description[sour[0]->type].parameter[i_item].name, 0)); value0 = *((long*)(sour[i_elem]->p_elem0+entity_description[sour[0]->type].parameter[i_item].offset)); sprintf(s, "%s0", entity_description[sour[0]->type].parameter[i_item].name); rpn_store(value0, NULL, rpn_create_mem(s, 0)); break; default: break; } #if DEBUG fprintf(stdout, " asserting value %e for %s\n", value, entity_description[sour[0]->type].parameter[i_item].name); fprintf(stdout, " asserting value %e for %s0\n", value0, entity_description[sour[0]->type].parameter[i_item].name); fflush(stdout); #endif } p_elem = targ[i_elem]->p_elem; rpn_clear(); /* push original value onto stack */ if (verbosity>1) fprintf(stdout, "prior value is %.15g for %s#%ld.%s at z=%.15gm\n", links->baseline_value[i_link][i_elem], links->target_name[i_link], targ[i_elem]->occurence, links->item[i_link], targ[i_elem]->end_pos); push_num(links->baseline_value[i_link][i_elem]); value = rpn(links->equation[i_link]); if (rpn_check_error()) exitElegant(1); rpn_clear(); if (value>links->maximum[i_link]) value = links->maximum[i_link]; if (value<links->minimum[i_link]) value = links->minimum[i_link]; if (verbosity>0) fprintf(stdout, "asserting value %.15g for %s#%ld.%s at z=%.15gm\n", value, links->target_name[i_link], targ[i_elem]->occurence, links->item[i_link], targ[i_elem]->end_pos); fflush(stdout); if (entity_description[elem_type].flags&HAS_LENGTH && entity_description[elem_type].parameter[param].offset==0) lengthChanged = 1; switch (data_type) { case IS_DOUBLE: *((double*)(p_elem+entity_description[elem_type].parameter[param].offset)) = value; break; case IS_LONG: *((long*)(p_elem+entity_description[elem_type].parameter[param].offset)) = nearestInteger(value); break; case IS_STRING: default: bombElegant("unknown/invalid variable quantity (assert_element_links)", NULL); exitElegant(1); } if (flags&LINK_ELEMENT_DEFINITION) change_defined_parameter_values(&targ[i_elem]->name, ¶m, &targ[i_elem]->type, &value, 1); if ((entity_description[targ[0]->type].parameter[param].flags&PARAM_CHANGES_MATRIX) && targ[i_elem]->matrix) { free_matrices(targ[i_elem]->matrix); tfree(targ[i_elem]->matrix); targ[i_elem]->matrix = NULL; compute_matrix(targ[i_elem], run, NULL); matrices_changed++; } } } #if DEBUG print_line(stdout, beamline); #endif if (lengthChanged) compute_end_positions(beamline); log_exit("assert_element_links"); return(matrices_changed); }
void tilt_matrices(VMATRIX *M, double tilt) { static VMATRIX Rot, IRot; static long initialized=0; VMATRIX Mr; double sin_tilt, cos_tilt; log_entry("tilt_matrices"); if (tilt==0) { log_exit("tilt_matrices"); return; } if (!initialized) { initialized = 1; initialize_matrices(&Rot, 1); initialize_matrices(&IRot, 1); } initialize_matrices(&Mr, M->order); if (fabs(tilt-PI/2)<1e-12) { sin_tilt = 1; cos_tilt = 0; } else if (fabs(tilt+PI/2)<1e-12) { sin_tilt = -1; cos_tilt = 0; } else if (fabs(tilt-PI)<1e-12 || fabs(tilt+PI)<1e-12) { sin_tilt = 0; cos_tilt = -1; } else { sin_tilt = sin(tilt); cos_tilt = cos(tilt); } /* Rotation matrix for (x, y) */ IRot.R[0][0] = Rot.R[0][0] = cos_tilt ; IRot.R[0][2] = -(Rot.R[0][2] = sin_tilt); IRot.R[2][0] = -(Rot.R[2][0] = -sin_tilt); IRot.R[2][2] = Rot.R[2][2] = cos_tilt ; /* Rotation matrix for (x', y') */ IRot.R[1][1] = Rot.R[1][1] = cos_tilt ; IRot.R[1][3] = -(Rot.R[1][3] = sin_tilt); IRot.R[3][1] = -(Rot.R[3][1] = -sin_tilt); IRot.R[3][3] = Rot.R[3][3] = cos_tilt ; IRot.R[4][4] = IRot.R[5][5] = Rot.R[4][4] = Rot.R[5][5] = 1; concat_matrices(&Mr, M, &Rot, 0); concat_matrices(M , &IRot, &Mr , 0); free_matrices(&Mr); log_exit("tilt_matrices"); }