Skip to content
147 changes: 147 additions & 0 deletions atintegrators/FeedbackPass.c
Original file line number Diff line number Diff line change
@@ -0,0 +1,147 @@
/* IdentityPass.c
Accelerator Toolbox
Revision 7/16/03
A.Terebilo
*/

#include "atelem.c"
#include "atlalib.c"
#ifdef MPI
#include <mpi.h>
#include <mpi4py/mpi4py.h>
#endif

struct elem
{
double GX;
double GY;
double GZ;
double *closed_orbit;
};

void FeedbackPass(double *r_in, int num_particles, struct elem *Elem)
{
int c;
double mx = 0.0;
double my = 0.0;
double mz = 0.0;
long npart = 0.0;

double *closed_orbit;
closed_orbit = Elem->closed_orbit;

double gx = Elem->GX;
double gy = Elem->GY;
double gz = Elem->GZ;

for (c = 0; c<num_particles; c++) { /*Loop over particles */
double *r6 = r_in+c*6;

if (!atIsNaN(r6[0])) {
npart += 1;
mx += r6[0];
my += r6[2];
mz += r6[5];
}
}

#ifdef MPI
MPI_Allreduce(MPI_IN_PLACE, &npart, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
MPI_Allreduce(MPI_IN_PLACE, &mx, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
MPI_Allreduce(MPI_IN_PLACE, &my, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
MPI_Allreduce(MPI_IN_PLACE, &mz, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
MPI_Barrier(MPI_COMM_WORLD);
#endif


if (npart>0.0){
mx /= npart;
my /= npart;
mz /= npart;
for (c = 0; c<num_particles; c++) { /*Loop over particles */
double *r6 = r_in+c*6;
if (!atIsNaN(r6[0])) {
r6[0] -= (mx-closed_orbit[0])*gx;
r6[2] -= (my-closed_orbit[2])*gy;
r6[5] -= (mz-closed_orbit[5])*gz;
}
}
}
}

#if defined(MATLAB_MEX_FILE) || defined(PYAT)
ExportMode struct elem *trackFunction(const atElem *ElemData,struct elem *Elem,
double *r_in, int num_particles, struct parameters *Param)
{
if (!Elem) {
double GX, GY, GZ, *closed_orbit;
GX=atGetOptionalDouble(ElemData,"Gx",0.0); check_error();
GY=atGetOptionalDouble(ElemData,"Gy",0.0); check_error();
GZ=atGetOptionalDouble(ElemData,"Gz",0.0); check_error();
closed_orbit = atGetDoubleArray(ElemData,"closed_orbit"); check_error();

Elem = (struct elem*)atMalloc(sizeof(struct elem));
Elem->GX=GX;
Elem->GY=GY;
Elem->GZ=GZ;
Elem->closed_orbit = closed_orbit;
}
FeedbackPass(r_in,num_particles,Elem);
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I think in other element we pass num_particles as the last argument but this is cosmetic

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

What I put here follows what is in BeamMomentsPass and SliceMomentsPass.

return Elem;
}

MODULE_DEF(FeedbackPass) /* Dummy module initialisation */

#endif /*defined(MATLAB_MEX_FILE) || defined(PYAT)*/


#ifdef MATLAB_MEX_FILE

void mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[])
{
if (nrhs >= 2) {

double *r_in;
const mxArray *ElemData = prhs[0];
int num_particles = mxGetN(prhs[1]);
struct elem El, *Elem=&El;

double GX,GY,GZ,*closed_orbit;

GX=atGetOptionalDouble(ElemData,"Gx",0.0); check_error();
GY=atGetOptionalDouble(ElemData,"Gy",0.0); check_error();
GZ=atGetOptionalDouble(ElemData,"Gz",0.0); check_error();
closed_orbit = atGetDoubleArray(ElemData,"closed_orbit");check_error();


Elem = (struct elem*)atMalloc(sizeof(struct elem));
Elem->GX=GX;
Elem->GY=GY;
Elem->GZ=GZ;
Elem->closed_orbit = closed_orbit;

if (mxGetM(prhs[1]) != 6) mexErrMsgIdAndTxt("AT:WrongArg","Second argument must be a 6 x N matrix: particle array");

/* ALLOCATE memory for the output array of the same size as the input */
plhs[0] = mxDuplicateArray(prhs[1]);
r_in = mxGetDoubles(plhs[0]);
FeedbackPass(r_in,num_particles,Elem);
}
else if (nrhs == 0) {
/* list of required fields */
plhs[0] = mxCreateCellMatrix(1,1);
mxSetCell(plhs[0],0,mxCreateString("closed_orbit"));
if(nlhs>1) /* optional fields */
{
plhs[1] = mxCreateCellMatrix(3,1);
mxSetCell(plhs[1],0,mxCreateString("Gx"));
mxSetCell(plhs[1],1,mxCreateString("Gy"));
mxSetCell(plhs[1],2,mxCreateString("Gz"));
}
}
else {
mexErrMsgIdAndTxt("AT:WrongArg","Needs 2 or 0 arguments");
}
}
#endif

37 changes: 36 additions & 1 deletion pyat/at/lattice/elements/basic_elements.py
Original file line number Diff line number Diff line change
Expand Up @@ -18,11 +18,12 @@
"SimpleRadiation",
"QuantumDiffusion",
"EnergyLoss",
"Feedback",
]

import warnings
from collections.abc import Iterable
from typing import Optional
from typing import Optional, Sequence

import numpy as np

Expand Down Expand Up @@ -627,4 +628,38 @@ def __init__(self, family_name: str, energy_loss: float, **kwargs):
super().__init__(family_name, EnergyLoss=energy_loss, **kwargs)


class Feedback(Element):
"""Transverse and Longitudinal Feedback element"""

def __init__(
self,
family_name: str,
Gx: float = 0.0,
Gy: float = 0.0,
Gz: float = 0.0,
closed_orbit: Sequence[float] = np.zeros(6),
**kwargs
):
"""
Args:
family_name: Name of the element
Gx: Feedback Gain in Horizontal (no units:
damping_time [turns] = 2 / Gx )
Gy: Feedback Gain in Vertical (no units:
damping_time [turns] = 2 / Gy )
Gz: Feedback Gain in Longitudinal (no units:
damping_time [turns] = 2 / Gz )
closed_orbit: 6D closed orbit to feedback onto

Default PassMethod: ``FeedbackPass``

Note that this element does not SET the closed orbit. That
is handled by the lattice (either full ring or linear maps for x and y,
or the ct for the longitudinal plane). Aan accurate closed orbit must be
provided in order to have a well behaving feedback.
"""
kwargs.setdefault("PassMethod", "FeedbackPass")
super().__init__(family_name, Gx=Gx, Gy=Gy, Gz=Gz, closed_orbit=closed_orbit, **kwargs)


Radiative.register(EnergyLoss)
10 changes: 7 additions & 3 deletions pyat/at/physics/fastring.py
Original file line number Diff line number Diff line change
Expand Up @@ -141,7 +141,9 @@ def simple_ring(
U0: float = 0.0,
name: str = "",
particle: str | Particle = "relativistic",
TimeLag: float | Sequence[float] = 0.0
TimeLag: float | Sequence[float] = 0.0,
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Isn't is possible to generate a simple ring from just an AT lattice?

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Well, that is typically what fast_ring is for. Simple_ring is more for inputs by hand. But you are right i could pass a lattice and make an even simpler fast_ring from a lattice. That is separate to this though.

orb6: Sequence[float] = np.zeros(6),

) -> Lattice:
"""Generates a "simple ring" based on a given dictionary
of global parameters
Expand Down Expand Up @@ -201,7 +203,8 @@ def simple_ring(
or a Particle object
TimeLag: Set the timelag of the cavities, Default=0. Can be scalar
or sequence of scalars (as with harmonic_number and Vrf).

orb6: 6d closed orbit. Needed for inclusion of feedbacks.

If the given emitx, emity or espread is 0, then no equlibrium emittance
is applied in this plane.
If the given tau is 0, then no radiation damping is applied for this plane.
Expand Down Expand Up @@ -265,7 +268,8 @@ def simple_ring(
# generate the linear tracking element, we set a length
# which is needed to give the lattice object the correct length
# (although it is not used for anything else)
lin_elem = M66("Linear", m66=Mat66, Length=circumference)

lin_elem = M66("Linear", m66=Mat66, Length=circumference, T1=-orb6, T2=orb6)

# Generate the simple radiation element
simplerad = SimpleRadiation(
Expand Down
Loading