/* * Note: The variables that save our state in this (and other fundtions) are * not local static variables, but are contained in the Lgm_MagModelInfo * structure that the user passes to it. This allows the function to be * thread-safe and re-entrant. */ double I_integrand( double s, _qpInfo *qpInfo ) { Lgm_Vector Bvec; int reset=1, done, Count; double f, g, Hremaining, Hdone, H, Htry, Hdid, Hnext, sgn=1.0, sdid, B, dS; Lgm_MagModelInfo *mInfo; /* * Get pointer to our auxilliary data structure. */ mInfo = (Lgm_MagModelInfo *)qpInfo; mInfo->Lgm_I_integrand_u_scale.x = 10.0; mInfo->Lgm_I_integrand_u_scale.y = 1.0; mInfo->Lgm_I_integrand_u_scale.z = 10.0; if ( mInfo->Lgm_I_integrand_JumpMethod == LGM_RELATIVE_JUMP_METHOD ) { /* * Set starting point. If this is the first call for this integral, * set point to the lower limit. */ if ( mInfo->Lgm_I_integrand_FirstCall == TRUE ) { mInfo->Lgm_I_integrand_FirstCall = FALSE; mInfo->Lgm_I_integrand_P = mInfo->Pm_South; mInfo->Lgm_I_integrand_S = 0.0; dS = s; } else { dS = s - mInfo->Lgm_I_integrand_S; } if (dS < 0.0) { H = -dS; sgn = -1.0; // H is a positive qnty } else { H = dS; sgn = 1.0; } } else if ( mInfo->Lgm_I_integrand_JumpMethod == LGM_ABSOLUTE_JUMP_METHOD ) { /* * This strategy starts at start each time. Slower(?), but seems to * limit roundoff error from accumulating? The speed increase of the * relative method really depends on how far apart the s points are * that the integrator is picking. If they are bouncing all over the * place the relative method still does alot of tracing. */ if ( mInfo->Lgm_I_integrand_FirstCall == TRUE ) { mInfo->Lgm_I_integrand_FirstCall = FALSE; } mInfo->Lgm_I_integrand_P = mInfo->Pm_South; mInfo->Lgm_I_integrand_S = 0.0; H = s; sgn = 1.0; // H is a positive qnty } else { printf("I_integrand: Error, unknown Jump Method ( Lgm_I_integrand_JumpMethod = %d\n", mInfo->Lgm_I_integrand_JumpMethod ); exit(-1); } /* * Get B-field at the given value of s. Need to advance along field line * by an amount H. May need to do more than one call to get there... * Setting Htry to be H is an attempt to make the full jump in one call to * MagStep. For very precise calculations, the number of jumps we do here * seems to be fairly critical. Making two jumps (i.e. Htry = H/2) and * limiting Hmax to 0.5 Re seems to work pretty well (perhaps its because * symetry between mirror points?). Could still play with the max step * size of 0.5 -- maybe something a bit higher would work just as well and * be more efficient? */ done = FALSE; Count = 0; Htry = 0.5*H; Hdone = 0.0; reset = TRUE; if (Htry > 0.5) Htry = 0.5; if ( fabs(mInfo->Lgm_I_integrand_S-s) < 1e-12 ) done = TRUE; while ( !done ) { Lgm_MagStep( &mInfo->Lgm_I_integrand_P, &mInfo->Lgm_I_integrand_u_scale, Htry, &Hdid, &Hnext, sgn, &sdid, &reset, mInfo->Bfield, mInfo ); Hdone += Hdid; // positive qnty mInfo->Lgm_I_integrand_S += sgn*Hdid; Hremaining = H - Hdone; if ( Count > 1000 ) { printf("I_integrand: Warning, early return because Count > 1000. Ill-conditioned integrand?\n"); done = TRUE; } else if ( fabs(mInfo->Lgm_I_integrand_S-s) < 1e-12 ){ done = TRUE; } else { if ( Htry > Hremaining ) Htry = Hremaining; } ++Count; } mInfo->Bfield( &mInfo->Lgm_I_integrand_P, &Bvec, mInfo ); B = Lgm_Magnitude( &Bvec ); /* * Compute integrand, ( 1 - B/Bm ) ^ (1/2) Make sure 1-B/Bm is not * negative. If it is, assume its zero. (Note that sometimes (due to * round-off error) it could be very slightly negative). */ g = 1.0 - B/mInfo->Bm; f = (g > 0.0) ? sqrt( g ) : 0.0; ++mInfo->Lgm_n_I_integrand_Calls; return( f ); }
int Lgm_TraceToSMEquat( Lgm_Vector *u, Lgm_Vector *v, double tol, Lgm_MagModelInfo *Info ) { Lgm_Vector u_scale; double Htry, Hdid, Hnext, Hmin, Hmax, s, sgn, r2; double Sa=0.0, Sb=0.0, Sc=0.0, d1, d2; double Za, Zb, Zc, Z; Lgm_Vector Ztmp; Lgm_Vector Pa, Pb, Pc, P; int done, reset; Hmax = 20.0; Hmin = 0.001; u_scale.x = 100.0; u_scale.y = 100.0; u_scale.z = 100.0; Za = Zb = Zc = 0.0; /* * Bracket the minimum. We want to find two points along * the field line such that the location of equatorial plane is * gauranteed to lie between them. To do this, we need to find * three points; Pa, Pb, and Pc such that; * * Pc > Pb > Pa * * and |z_sm( Pa )| > |z_sm( Pb )| * and |z_sm( Pc )| > |z_sm( Pb )| * * */ /* * Set the start point, Pa and Za (the absolute value of the SM z coord). */ Pa = *u; Lgm_Convert_Coords( &Pa, &Ztmp, GSM_TO_SM, Info->c ); Za = fabs( Ztmp.z ); Sa = 0.0; /* * If we are above the SM plane (i.e. in the northern hemisphere), lets * try tracing against the field direction. */ sgn = ( Ztmp.z > 0.0 ) ? -1.0 : 1.0; /* * Keep stepping along the field line until we have a |Zsm| less than the start point. */ P = Pa; done = FALSE; Htry = Za; while ( !done ) { //if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, 1.0e-7, sgn, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, sgn, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); Lgm_Convert_Coords( &P, &Ztmp, GSM_TO_SM, Info->c ); Z = fabs( Ztmp.z ); if ( (P.x > Info->OpenLimit_xMax) || (P.x < Info->OpenLimit_xMin) || (P.y > Info->OpenLimit_yMax) || (P.y < Info->OpenLimit_yMin) || (P.z > Info->OpenLimit_zMax) || (P.z < Info->OpenLimit_zMin) ) { /* * Open FL! */ v->x = v->y = v->z = 0.0; return(0); } else if ( Z < Za ) { done = TRUE; Pb = P; Zb = Z; Sb += Hdid; } else { Pa = P; Za = Z; Sa = 0.0; } Htry = Z; } /* * Keep stepping along the field line until we complete the bracket triple. */ done = FALSE; while (!done) { P = Pb; //if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, 1.0e-7, sgn, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, sgn, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); Lgm_Convert_Coords( &P, &Ztmp, GSM_TO_SM, Info->c ); Z = fabs( Ztmp.z ); if ( (P.x > Info->OpenLimit_xMax) || (P.x < Info->OpenLimit_xMin) || (P.y > Info->OpenLimit_yMax) || (P.y < Info->OpenLimit_yMin) || (P.z > Info->OpenLimit_zMax) || (P.z < Info->OpenLimit_zMin) ) { /* * Open FL! */ v->x = v->y = v->z = 0.0; return(0); } else if ( Z < Zb ) { Pa = Pb; Za = Zb; Sa = 0.0; Pb = P; Zb = Z; Sb = Hdid; r2 = P.x*P.x + P.y*P.y + P.z*P.z; Htry = (Hnext < Hmax) ? ((Hnext > Hmin) ? Hnext : Hmin) : Hmax; if (Htry*Htry > r2) Htry = 0.25*sqrt(r2); } else { Pc = P; Zc = Z; Sc = Sb + Hdid; done = TRUE; } } /* * We have a bracket. Now go in for the kill. * Use golden-section search to converge toward minimum. * (Sa, Sb, Sc) are the distances of the triple points along * the FL. For convenience, we maintain Sa = 0.0, so that * Sb and Sc are always the distances relative to Pa. (And * so Sc is always the width of the bracketed interval. So when * Sc gets small enough we bail out and take Pb as the min). */ done = FALSE; while (!done) { d1 = Sb; d2 = Sc - Sb; if ( Sc < tol ) { done = TRUE; } else if ( d1 > d2 ) { P = Pa; Htry = 0.381966011*d1; //if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, 1.0e-7, sgn, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, sgn, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); Lgm_Convert_Coords( &P, &Ztmp, GSM_TO_SM, Info->c ); Z = fabs( Ztmp.z ); if ( Z < Zb ) { Pb = P; Zb = Z; Sb = Hdid; } else { Pa = P; Za = Z; Sb -= Hdid; Sc -= Hdid; } } else { P = Pb; Htry = 0.381966011*d2; //if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, 1.0e-7, sgn, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, sgn, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); Lgm_Convert_Coords( &P, &Ztmp, GSM_TO_SM, Info->c ); Z = fabs( Ztmp.z ); if ( Z < Zb ) { Pb = P; Zb = Z; Sb += Hdid; } else { Pc = P; Zc = Z; Sc = Sb + Hdid; } } } *v = Pb; return( 1 ); }
int Lgm_TraceToYZPlane( Lgm_Vector *u, Lgm_Vector *v, double Xtarget, double sgn0, double tol, Lgm_MagModelInfo *Info ) { Lgm_Vector u_scale; double Htry, Hdid, Hnext, Hmin, Hmax, s, sgn, fsgn0; double Sa, Sb, B, f, r2, z2, r3, L, R, hhh; double Stotal; Lgm_Vector Btmp; Lgm_Vector Pa, Pb, P; int done, reset=TRUE, bracketed=FALSE; long int Ntotal; //printf("\n\n\n***************************************\n"); Pb.x = Pb.y = Pb.z = 0.0; Hmax = 20.0; Hmin = 0.01; u_scale.x = 100.0; u_scale.y = 100.0; u_scale.z = 100.0; /* * Set the start point, Pa */ Pa = *u; f = Xtarget - Pa.x; fsgn0 = ( f<0.0) ? -1.0 : 1.0; Sa = 0.0; Stotal = 0.0; Ntotal = 0; /* * Choose a good step size to try. */ Lgm_Convert_Coords( &Pa, &P, GSM_TO_SM, Info->c ); R = Lgm_Magnitude( &Pa ); r2 = R*R; z2 = P.z*P.z; L = 9e99; if ( fabs( r2 - z2 ) < 1e-2 ) { /* * High L */ Htry = 2.0; } else { r3 = r2*R; L = r3 / (r2 - z2); Htry = ( L < 4.5 ) ? 1.165*L : 4.5; } if ( Htry < 1e-4 ) Htry = 0.001; //Htry = 0.01; /* * Now, bracket minimum. And do a bisection search for the minimum. */ done = FALSE; sgn = sgn0; P = Pa; Sb = -9e99; //printf("Htry =%g sgn = %g\n", Htry, sgn ); //printf("Xtarget = %g f = %g P = %g %g %g\n", Xtarget, f, P.x, P.y, P.z); while ( !done ) { //if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, 1.0e-7, sgn, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, sgn, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); Stotal += Hdid; ++Ntotal; R = Lgm_Magnitude( &P ); r2 = R*R; r3 = r2*R; L = r3 / (r2 - z2); f = Xtarget - P.x; //printf("Xtarget = %g f = %g P = %g %g %g\n", Xtarget, f, P.x, P.y, P.z); if ( (P.x > Info->OpenLimit_xMax) || (P.x < Info->OpenLimit_xMin) || (P.y > Info->OpenLimit_yMax) || (P.y < Info->OpenLimit_yMin) || (P.z > Info->OpenLimit_zMax) || (P.z < Info->OpenLimit_zMin) || ( s > 1000.0 ) ) { /* * Open FL! */ //v->x = v->y = v->z = 0.0; *v = P; //printf("OPEN!\n\n\n"); return(0); } else if ( r2 < 1.0 ) { return( -1 ); } else if ( (Stotal > 300.0) || ( Ntotal > 1500 ) ) { printf("Lgm_TraceToYZPlane: %s:%d Field line too long or too many iterations: Stotal / Ntotal: %g / %ld\n", __FILE__, __LINE__, Stotal, Ntotal ); return( 0 ); } else if ( fabs( f ) < tol ){ done = TRUE; } else if ( fsgn0*f > 0.0 ){ Pa = P; Sa += Hdid; sgn = sgn0; if ( bracketed ) { Htry = fabs( Sb - Sa )/2.0; } else if ( L < 4.5 ) { hhh = 1.165*L+1.0 - R + 1.0; //if ( hhh > 1e-4 ) Htry = hhh; if ( hhh > R-1.0 ) { hhh = (R-1.0)/2.0; } } /* printf("A: Sa, Sb, P = %g %g (%g, %g, %g) Htry = %f\n", Sa, Sb, P.x, P.y, P.z, Htry); */ } else { Pb = P; if ( bracketed ) Sb += sgn*sgn0*Hdid; else Sb = Sa + Hdid; bracketed = TRUE; sgn = -sgn0; Htry = fabs( Sb - Sa )/2.0; /* printf("B: Sa, Sb, P = %g %g (%g, %g, %g) Htry = %f\n", Sa, Sb, P.x, P.y, P.z, Htry); */ } } /* * */ *v = Pb; //printf("B: Sa, Sb, P = %g %g (%g, %g, %g) Htry = %f\n", Sa, Sb, P.x, P.y, P.z, Htry); // if ( fabs( Xtarget - Pb.x ) > 2.0*tol ) return( -1 ); return( 1 ); }
int Lgm_TraceToSphericalEarth( Lgm_Vector *u, Lgm_Vector *v, double TargetHeight, double sgn, double tol, Lgm_MagModelInfo *Info ) { Lgm_Vector u_scale; double Htry, Htry_max, Hdid, Hnext, Hmin, Hmax, s; double Sa=0.0, Sc=0.0, d; double Rinitial, Fa, Fb, Fc, F; double Height, StartHeight; double Height_a, Height_b, Height_c, HeightPlus, HeightMinus, direction; Lgm_Vector Pa, Pc, P; int done, reset, AboveTargetHeight, Count; reset = TRUE; Info->Trace_s = 0.0; Sa = Sc = 0.0; /* * Determine our initial geocentric radius in km. (u is assumed to be in * units of Re where we define Re to be WGS84_A.) */ Rinitial = WGS84_A*Lgm_Magnitude( u ); // km /* * The Earth is a spheroid. It is more flattened at the poles than at the * equator. Check to see if we are initially at a height where we need to * worry about it. For WGS84 Earth shape model, the equatorial radius is * WGS84_A and polar radius is WGS84_B (WGS84_B is about 20km smaller than * WGS84_A). * * Note that although we are trying to trace to a spherical Earth in this * routine, we still dont want to drop below the surface of the real Earth * (IGRF doesnt like that so much). * */ if ( Rinitial < WGS84_B ) { // must be inside the Earth, which is no good -- bail with // LGM_INSIDE_EARTH error code return( LGM_INSIDE_EARTH ); } else { // We are at least at or above WGS84_B. We could still be in trouble in // terms of being inside the Earth, so we have to check more // thouroughly now. Determine Geodetic Height Lgm_WGS84_to_GeodHeight( u, &Height ); if ( Height < 0.0 ) { // inside the Earth, which is no good -- bail with error return( LGM_INSIDE_EARTH ); } } /* * If we get here, the initial point must be above the surface of the * (spheroidal) Earth. * * Now check to see if we are currently above or below the target height. * (reset Height to be distance in km above the spherical approx of the Earth.) */ Height = Rinitial - WGS84_A; // distance above spherical Earth AboveTargetHeight = ( Height > TargetHeight ) ? TRUE : FALSE; /* * Initialize some things */ Pa.x = Pa.y = Pa.z = 0.0; Pc.x = Pc.y = Pc.z = 0.0; P = *u; Hmax = Info->Hmax; Hmax = 1.0; Hmin = 0.001; Hmin = 1e-8; u_scale.x = u_scale.y = u_scale.z = 1.0; Height = Height_a = Height_b = Height_c = 0.0; F = Fa = Fb = Fc = 0.0; //printf("\nHmax = %g\n", Hmax); /* * Save initial point if we need to */ if (Info->SavePoints) fprintf(Info->fp, "%f \t%f\t %f\t 3\n", u->x, u->y, u->z); /* * If we are above the target height, a potential problem will occur if the * field line is open in the direction of tracing. Fortunately, that * situation is easy to detect and is taken care of later... * * Here we need to test for and deal with the opposite problem. If we are * already below the target height, it is possible that the field line * never gets high enough to reach the target height. This is more tricky * to deal with because we can run the risk of getting below the surface of * the Earth. We need to try to get back above the target height to test * for this. If we cannot get above we need to bail out of the routine * altogether. */ StartHeight = Height; if ( !AboveTargetHeight ) { /* * Determine which direction along the field line will take us higher. */ Htry = 0.001; // sgn = +1 P = *u; if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, 1.0, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); HeightPlus = WGS84_A*(Lgm_Magnitude( &P )-1.0); // sgn = -1 P = *u; if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, -1.0, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); HeightMinus = WGS84_A*(Lgm_Magnitude( &P )-1.0); direction = ( HeightPlus > HeightMinus ) ? 1.0 : -1.0; /* * We are already at or below target height. Trace until we are not. */ done = FALSE; //reset = TRUE; while ( !done ) { Htry = fabs(0.9*(TargetHeight - Height)); // This computes Htry as 90% of the distance to the TargetHeight if (Htry > 0.1) Htry = 0.1; // If its bigger than 0.1 reset it to 0.1 -- to be safe. if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, direction, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); //printf("1. PPPPPPPPPP = %g %g %g\n", P.x, P.y, P.z); Sa += Hdid; Info->Trace_s += Hdid; Height = WGS84_A*(Lgm_Magnitude( &P )-1.0); F = Height - TargetHeight; if ( F > 0.0 ) { done = TRUE; } else if ( Height < StartHeight ) { // We are going back down again -- Target Height unreachable? -- Bail out return( LGM_TARGET_HEIGHT_UNREACHABLE ); } } } /* * Bracket the zero first. * We want to stop in the ionosphere at a certain height above the Earth * (assumed to be a spehere of radius WGS84_A in this routine). We need * to find 2 points along the field line such that the intersection of the * FL with the sphere is guaranteed to lie between them. To do this, we * need to find two points; Pa and Pc such that; * * and Height( Pa ) - TargetHeight > 0.0 * and Height( Pc ) - TargetHeight < 0.0 * * I.e., F = Height - TargetHeight has opposite signs * * Set the start point, Pa and Fa. (Fa is difference between Height_a and * TargetHeight.) * */ Pa = P; Height_a = WGS84_A*(Lgm_Magnitude( &Pa )-1.0); Fa = Height_a - TargetHeight; /* * Get an initial Htry that is safe -- i.e. start off slowly * We dont really know where we are, so be conservative on the first try. * If if ( Lgm_MagStep() gives back an Hnext thats higher, we'll crank Htry up then... */ Htry = 0.9*Height_a; // This computes Htry as 90% of the distance to the Earth's surface (could be small if we are already close!) if (Htry > Hmax) Htry = Hmax; // If its bigger than Hmax reset it to Hmax -- to be safe. /* * Keep stepping along the field line until we drop below the target * height, TargetHeight. (Or bail if its open). This completes the endpoints of the * bracket pair. We get these points first, because there are field lines * that have more than one local minimum in fabs(Height-TargetHeight). So find Pa and Pc * first and then complete the triple later. CHECK on this. I changed this * from a minima search to a simple root finder. Is it OK still? */ P = Pa; done = FALSE; //reset = TRUE; Count = 0; while ( !done ) { if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, sgn, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); Height = WGS84_A*(Lgm_Magnitude( &P )-1.0); //Lgm_Vector BBBBB; //Info->Bfield( &P, &BBBBB, Info ); //printf("2. PPPPPPPPPP = %g %g %g BBBBBBBBBB = %g %g %g Height = %g HTRY = %g HDID = %g HNEXT = %g\n", P.x, P.y, P.z, BBBBB.x, BBBBB.y, BBBBB.z, Height, Htry, Hdid, Hnext); F = Height - TargetHeight; if ((F > 0.0) && (Info->SavePoints)) fprintf(Info->fp, "%f \t%f\t %f\t 2\n", P.x, P.y, P.z); if ( (P.x > Info->OpenLimit_xMax) || (P.x < Info->OpenLimit_xMin) || (P.y > Info->OpenLimit_yMax) || (P.y < Info->OpenLimit_yMin) || (P.z > Info->OpenLimit_zMax) || (P.z < Info->OpenLimit_zMin) || ( s > 1000.0 ) ) { /* * Open FL! */ v->x = v->y = v->z = 0.0; return(0); } else if ( F < 0.0 ) { done = TRUE; Pc = P; Fc = F; Height_c = Height; //Sc += Hdid; Sc = Sa + Hdid; } else { Pa = P; Fa = F; Height_a = Height; Sa += Hdid; } //if ( fabs(Hdid-Htry) > 1e-6 ) Htry = Hnext; // adaptively reset Htry Htry = Hnext; // adaptively reset Htry //printf("A. Htry = %g\n", Htry); /* * Go no farther than some small distance below * the target Height. Also respect Hmin and Hmax. */ Htry_max = 0.9*Height; if (Htry > Htry_max) Htry = Htry_max; //printf("B. Htry = %g Htry, Hmin, Hmax, Htry_max = %g %g %g %g\n", Htry, Htry, Hmin, Hmax, Htry_max); if (Htry < Hmin) Htry = Hmin; else if (Htry > Hmax) Htry = Hmax; //printf("C. Htry = %g\n", Htry); if ( Count > 1000) { printf("File: %s Lgm_TraceToSphericalEarth(), Line: %d; Too many iterations trying to reach target height (are we in a weird field region?) Returning with -1.\n", __FILE__, __LINE__ ); //exit(0); return(-1); } ++Count; } if ( ((Fc > 0.0) && (Fa > 0.0)) || ((Fc < 0.0) && (Fa < 0.0)) ) { // No bracket return(0); } /* * We have a bracket. Now go in for the kill. * * (Note: If all we wanted to know is whether or not the line hits the * TargetHeight, we could stop here: it must or we wouldnt have a minimum * bracketed.) * * Use golden-section search to converge toward minimum. * (Sa, Sb, Sc) are the distances of the triple points along * the FL. For convenience, we maintain Sa = 0.0, so that * Sb and Sc are always the distances relative to Pa. (And * so Sc is always the width of the bracketed interval. So when * Sc gets small enough we bail out and take Pb as the min). * */ //reset = TRUE; if (1==1) { done = FALSE; while (!done) { d = Sc - Sa; //if ( fabs(F) < tol ) { if ( fabs(d) < tol ) { done = TRUE; } else { P = Pa; //Htry = LGM_1M_1O_GOLD*d; // LGM_1M_1O_GOLD is 0.381966... // if ( Htry > 2.0*fabs(F) ) { //Htry = LGM_1M_1O_GOLD*d; // LGM_1M_1O_GOLD is 0.381966... Htry = 0.5*fabs(d); // LGM_1M_1O_GOLD is 0.381966... // } if ( Lgm_MagStep( &P, &u_scale, Htry, &Hdid, &Hnext, sgn, &s, &reset, Info->Bfield, Info ) < 0 ) return(-1); //printf("3. PPPPPPPPPP = %g %g %g HTRY = %g\n", P.x, P.y, P.z, Htry); Height = WGS84_A*(Lgm_Magnitude( &P )-1.0); F = Height - TargetHeight; if ( F >= 0.0 ) { Pa = P; Fa = F; Sa += Hdid; } else { Pc = P; Fc = F; Sc = Sa + Hdid; } } } Info->Trace_s = Sa; /* * Take average as the final answer. */ v->x = 0.5*(Pa.x + Pc.x); v->y = 0.5*(Pa.y + Pc.y); v->z = 0.5*(Pa.z + Pc.z); //printf( "Lgm_TraceToSphericalEarth: v = %g %g %g\n", v->x, v->y, v->z); } if (0==1) { /* * Try Brent's method */ //printf("Sa, Sb = %g %g Fa, Fb = %g %g tol = %g\n", Sa, Sb, Fa, Fb, tol); double Sz, Fz; Lgm_Vector Pz; BrentFuncInfoP f; f.u_scale = u_scale; f.Htry = Htry; f.sgn = sgn; f.reset = reset; f.Info = Info; f.func = &seFunc; f.Val = TargetHeight; Lgm_zBrentP( Sa, Sc, Fa, Fc, Pa, Pc, &f, tol, &Sz, &Fz, &Pz ); Fc = Fz; Sc = Sz; Pc = Pz; //printf("Sa, Sb = %g %g Fa, Fb = %g %g tol = %g\n", Sa, Sb, Fa, Fb, tol); v->x = Pz.x; v->y = Pz.y; v->z = Pz.z; Info->Trace_s = Sz; } if (Info->SavePoints) fprintf(Info->fp, "%f \t%f\t %f\t 2\n", v->x, v->y, v->z); if ( Info->VerbosityLevel > 2 ) printf("Lgm_TraceToSphericalEarth(): Number of Bfield evaluations = %ld\n", Info->Lgm_nMagEvals ); return( 1 ); }