// calcuates Dij vector: It points outer of domain
void DijVector(pFace face, pVertex &oppositeVertex, double *Dij){
	int i;
	double xyz[3][3], xyzOV[3];

	// get all vertices coordinates
	for (i=0 ;i<3; i++){
		V_coord( (pVertex)face->get(0,i), xyz[i] );
	}
	// get opposite vertex to boundary face coordinates
	V_coord( oppositeVertex, xyzOV );

	double a[3] = { xyz[1][0]-xyz[0][0], xyz[1][1]-xyz[0][1], xyz[1][2]-xyz[0][2] } ;	// vector over face's edge
	double b[3] = { xyz[2][0]-xyz[0][0], xyz[2][1]-xyz[0][1], xyz[2][2]-xyz[0][2] } ;	// vector over another face's edge
	double c[3] = { xyzOV[0]-xyz[0][0], xyzOV[1]-xyz[0][1], xyzOV[2]-xyz[0][2] } ;		// vector from face's center to opposite vertex

	// Dij: Dij = a x b
	double _Dij[3]={.0,.0,.0};
	computeCrossProduct(a,b,_Dij);

	// if inner product between Dij and c vector is negative, it means Dij points inside domain
	if ( computeDotProduct( c, _Dij ) > .0 ){
		for (i=0;i<3;i++){
			_Dij[i] = -_Dij[i];
		}
	}
	for (i=0;i<3;i++){
		Dij[i] = _Dij[i]/6.0;
	}
}
Exemplo n.º 2
0
//------------------------------------------------------------------------------
void Solver::computeAcceleration (unsigned char res)
{
    IDList::const_iterator i = mFluidParticles[res]->ActiveIDs.begin();
    IDList::const_iterator e = mFluidParticles[res]->ActiveIDs.end();

    for (; i != e; i++)
    {
        float *pos = &mFluidParticles[res]->Positions[2*(*i)];
        float *vel = &mVelocities[res][2*(*i)];
        float *acc = &mAccelerations[res][2*(*i)];
        float den = mDensities[res][*i];
        float pre = mPressures[res][*i];


        //	std::cout << mDensities[LOW][*i] << std::endl;

        float ene[2];
        ene[0] = 0.0f;
        ene[1] = 0.0f;
        float psiSum = 0.0f;

        // reset acc
        acc[0] = 0.0f;
        acc[1] = 0.0f;

        // reserve mem for tension force
        float accT[2];
        accT[0] = 0.0f;
        accT[1] = 0.0f;

        // reserve mem for boundary force
        float accB[2];
        accB[0] = 0.0f;
        accB[1] = 0.0f;

        //======================================================================
        // 	Compute force contribution from the same domain
        //======================================================================

        // get neighbor ids
        mFluidHashTable[res]->Query(Vector2f(pos));
        const IDList& neighbors = mFluidHashTable[res]->GetResult();

        IDList::const_iterator j = neighbors.begin();
        IDList::const_iterator je = neighbors.end();

        for (; j != je; j++)
        {
            float *posj = &mFluidParticles[res]->Positions[2*(*j)];
            float *velj = &mVelocities[res][2*(*j)];
            float denj = mDensities[res][*j];
            float prej = mPressures[res][*j];
            float dist = computeDistance(pos, posj);
            float xij[2];
            xij[0] = pos[0] - posj[0];
            xij[1] = pos[1] - posj[1];
            float vij[2];
            vij[0] = vel[0] - velj[0];
            vij[1] = vel[1] - velj[1];

            if (dist < mConfiguration.EffectiveRadius[res])
            {
                // compute SPH pressure coefficient
                float coeff = (pre/(den*den) + prej/(denj*denj));
                float vx = computeDotProduct(xij, vij);
                float h = mConfiguration.EffectiveRadius[res];

                // check if artificial viscosity shoud be applied
                if (vx < 0.0f)
                {
                    // add artificial velocity to the SPH Force coefficient
                    float x2 = dist*dist;
                    float nu = -2.0f*mConfiguration.Alpha*h*
                               mConfiguration.SpeedSound/(den + denj);
                    coeff += nu*vx/(x2 + 0.01f*h*h);
                }

                // evaluate kernel gradient and add contribution of particle j
                // to acc
                Vector2f grad = evaluateKernelGradient(Vector2f(xij), dist, h);
                acc[0] += mBlendValues[res][*j]*coeff*grad.X;
                acc[1] += mBlendValues[res][*j]*coeff*grad.Y;

                // compute tension force and add to the acc
                float kernelT = mBlendValues[res][*j]*evaluateKernel(dist, h);
                accT[0] -= mConfiguration.TensionCoefficient*xij[0]*kernelT;
                accT[1] -= mConfiguration.TensionCoefficient*xij[1]*kernelT;


                float w2 = mConfiguration.FluidParticleMass[res]/
                           mConfiguration.FluidParticleMass[LOW];
                float mw = w2*mexicanHat2D
                           (
                               xij[0]/mConfiguration.EffectiveRadius[LOW],
                               xij[1]/mConfiguration.EffectiveRadius[LOW]
                           );
                ene[0] += velj[0]*mw;
                ene[1] += velj[1]*mw;
                psiSum += mw;

                mStates[res][*i] |= (mStates[res][*j] << 1) & 0x04;

            }

        }

        //======================================================================
        // 	Compute force contribution from the contemp. domain
        //======================================================================

        float accC[2];
        accC[0] = 0.0f;
        accC[1] = 0.0f;

        float accCT[2];
        accCT[0] = 0.0f;
        accCT[1] = 0.0f;

        unsigned char resc = (res + 1) % 2;
        mFluidHashTable[resc]->Query
        (
            Vector2f(pos),
            mConfiguration.EffectiveRadius[LOW]
        );
        const IDList& neighborsc = mFluidHashTable[resc]->GetResult();

        IDList::const_iterator jc = neighborsc.begin();
        IDList::const_iterator jce = neighborsc.end();

        for (; jc != jce; jc++)
        {
            float *posj = &mFluidParticles[resc]->Positions[2*(*jc)];
            float *velj = &mVelocities[resc][2*(*jc)];
            float denj = mDensities[resc][*jc];
            float prej = mPressures[resc][*jc];
            float dist = computeDistance(pos, posj);
            float xij[2];
            xij[0] = pos[0] - posj[0];
            xij[1] = pos[1] - posj[1];
            float vij[2];
            vij[0] = vel[0] - velj[0];
            vij[1] = vel[1] - velj[1];

            if (dist < mConfiguration.EffectiveRadius[LOW])
            {
                // compute SPH pressure coefficient
                float coeff = (pre/(den*den) + prej/(denj*denj));
                float vx = computeDotProduct(xij, vij);
                float h = mConfiguration.EffectiveRadius[res];
                float hc = mConfiguration.EffectiveRadius[resc];

                if (vx < 0.0f)
                {
                    // add artificial velocity to the SPH Force coefficient
                    float x2 = dist*dist;
                    float nu = -2.0f*mConfiguration.Alpha*h*
                               mConfiguration.SpeedSound/(den + denj);
                    coeff += nu*vx/(x2 + 0.01f*h*h);
                }

                Vector2f grad = evaluateKernelGradient
                                (
                                    Vector2f(xij),
                                    dist,
                                    h
                                );
                Vector2f gradc = evaluateKernelGradient
                                 (
                                     Vector2f(xij),
                                     dist,
                                     hc
                                 );

                accC[0] += mBlendValues[resc][*jc]*coeff*
                           (grad.X + gradc.X)*0.5f;
                accC[1] += mBlendValues[resc][*jc]*coeff*
                           (grad.Y + gradc.Y)*0.5f;

                float wt = evaluateKernel(dist, h);
                float wct = evaluateKernel(dist, hc);

                accCT[0] -= mConfiguration.TensionCoefficient*xij[0]*
                            (wt + wct)/2.0f*mBlendValues[resc][*jc];
                accCT[1] -= mConfiguration.TensionCoefficient*xij[1]*
                            (wt + wct)/2.0f*mBlendValues[resc][*jc];


                float w3 = 1.0f;
                if (res == HIGH && dist > mConfiguration.EffectiveRadius[HIGH])
                {
                    w3 = 0.0f;
                }

                float w2 = mConfiguration.FluidParticleMass[res]/
                           mConfiguration.FluidParticleMass[LOW];
                float mw = w3*w2*mexicanHat2D
                           (
                               xij[0]/mConfiguration.EffectiveRadius[LOW],
                               xij[1]/mConfiguration.EffectiveRadius[LOW]
                           );
                ene[0] += velj[0]*mw;
                ene[1] += velj[1]*mw;
                psiSum += mw;

                mStates[res][*i] |= (mStates[res][*jc] << 1) & 0x04;
            }

        }

        //======================================================================
        // 	compute penalty force of the boundary
        //======================================================================

        mBoundaryHashTable->Query(Vector2f(pos));
        const IDList& bneighbors = mBoundaryHashTable->GetResult();

        IDList::const_iterator k = bneighbors.begin();
        IDList::const_iterator ke = bneighbors.end();

        for (; k != ke; k++)
        {
            float *posk = &mBoundaryParticles->Positions[2*(*k)];
            float xik[2];
            xik[0] = pos[0] - posk[0];
            xik[1] = pos[1] - posk[1];
            float dist = computeNorm(xik);

            if (dist < mConfiguration.EffectiveRadius[LOW])
            {
                float w = evaluateBoundaryWeight
                          (
                              dist,
                              mConfiguration.EffectiveRadius[LOW],
                              mConfiguration.SpeedSound
                          );
                float c = mConfiguration.BoundaryParticleMass/
                          (mConfiguration.FluidParticleMass[LOW] +
                           mConfiguration.BoundaryParticleMass);

                float d = c*w;
                accB[0] += d*xik[0];
                accB[1] += d*xik[1];
            }

        }

        //======================================================================
        // 	Update acceleration of particle i
        //======================================================================

        acc[0] *= -mConfiguration.FluidParticleMass[res];
        acc[1] *= -mConfiguration.FluidParticleMass[res];
        accC[0] *= -mConfiguration.FluidParticleMass[resc];
        accC[1] *= -mConfiguration.FluidParticleMass[resc];
        accCT[0] *= (mConfiguration.FluidParticleMass[resc]/
                     mConfiguration.FluidParticleMass[res]);
        accCT[1] *= (mConfiguration.FluidParticleMass[resc]/
                     mConfiguration.FluidParticleMass[res]);
        acc[0] += accT[0];
        acc[1] += accT[1];
        acc[0] += accC[0];
        acc[1] += accC[1];
        acc[0] += accCT[0];
        acc[1] += accCT[1];
        acc[0] += accB[0];
        acc[1] += accB[1];
        acc[1] -= 9.81f;

        float energy = 1.0f/(psiSum*psiSum*mConfiguration.EffectiveRadius[res])*
                       (ene[0]*ene[0] + ene[1]*ene[1]);


        if (pos[0] > 0.5f)
        {
            //mFluidParticles[res]->Colors[*i] = 1.0f;
            mStates[res][*i] |= 0x08;
        }
        else
        {
            //mFluidParticles[res]->Colors[*i] = 0.0f;
        }
        mFluidParticles[res]->Colors[*i] = std::min
                                           (
                                               abs(energy)/200.0f,
                                               1.0f
                                           );

//		xicm[0] = pos[0] - xicm[0]/mt;
//		xicm[1] = pos[1] - xicm[1]/mt;
//
//		float col = std::min
//		(
//			computeNorm(xicm), 0.003f
//		)/0.003f;
//
//		std::cout << col << std::endl;
//
//		mFluidParticles[res]->Colors[*i] = col;
    }
}
Exemplo n.º 3
0
void MlMaximumEntropyModel::learnLMBFGS(MlTrainingContainer* params)
{
	params->initialize();

	const MlDataSet* trainingDataSet = params->getTrainingDataSet();
	const MlDataSet* testDataSet     = params->getTestDataSet();
	const vector<size_t>& trainingIdxs = params->getTrainingIdxs();
	const vector<size_t>& testIdxs	   = params->getTestIdxs();
	const size_t numSamples = trainingIdxs.size();
	const bool performTest = (testDataSet && testIdxs.size()>0);

	if (trainingDataSet->getNumClasess() < 2)
		error("learnLMBFGS accepts only datasets with 2 or more classes, yor data has ",
			trainingDataSet->getNumClasess());

	const double  lambda  = params->getLambda();
	const size_t memorySize = params->getLmBfgsMemorySize(); 
	const size_t reportFrequency = params->getVerboseLevel();
	const double perplexityDelta  = params->getPerplexityDelta();
	const size_t numClasses  = trainingDataSet->getNumClasess();
	const size_t numFeatures = trainingDataSet->getNumBasicFeatures(); // F
	const size_t numTraining = trainingIdxs.size();    // N

	const size_t numRestarts=3;

	// data structures used for training
	vector< vector<double> >& w = weights_; // class X features
	vector< vector<double> > wOld(numClasses, vector<double>(numFeatures,0.0)); // class X features
	vector< vector<double> > wtx(numClasses, vector<double>(numTraining, 0.0));  // class X samples
	vector< vector<double> > qtx(numClasses, vector<double>(numTraining, 0.0));  // class X samples
	vector< vector<double> > q(numClasses, vector<double>(numFeatures,0.0));       // class X features
	vector< vector<double> > g(numClasses, vector<double>(numFeatures,0.0));       // class X features
	vector< vector<double> > gOld(numClasses, vector<double>(numFeatures,0.0));  // class X features

	vector< vector<float> > trainingProbs(numClasses, vector<float>(numTraining));
	vector< vector<float> > testProbs(numClasses, vector<float>(numTraining));
	vector< vector<double> > bestW(numClasses, vector<double>(numFeatures));

	// initialize weights
	if (params->getInputPath().length() > 1)
	{
		const string modelFile = params->getInputPath() + "_scr.txt";
		if (readModel(modelFile.c_str()))
			params->setIndClearWeights(false);
	}

	if (params->getIndClearWeights())
		weights_.clear();

	weights_.resize(numClasses, vector<double>(numFeatures,0.0));

	double previousPerplexity = MAX_FLOAT;
	float  bestTestError=1.0;
	size_t bestTestRound=0;
	float  bestTrainingError=1.0;
	size_t bestTrainingRound=0;

	bool terminateTraining = false;
	size_t totalRounds=0;
	size_t megaRound=0;
	for ( megaRound=0; megaRound<numRestarts; megaRound++)
	{
		// first round
		computeGradient(trainingDataSet, trainingIdxs, w, wtx, lambda, g);

		const double gtg = computeDotProduct(g,g);
		const double denominator = 1.0 / sqrt(gtg);
		for (size_t c=0; c<numClasses; c++)
			for (size_t i=0; i<numFeatures; i++)
				q[c][i]=g[c][i]*denominator;

		// qtx <- qTx
		for (size_t c=0; c<numClasses; c++)
			for (size_t i=0; i<numSamples; i++)
			{
				const MlSample& sample = trainingDataSet->getSample(trainingIdxs[i]);
				qtx[c][i]=computeDotProduct(q[c],sample.pairs);
			}

		// eta <- lineSearch(...)
		double eta = lineSearch(trainingDataSet, trainingIdxs, w, wtx, qtx, g, q, lambda);
		//cout << "eta = " << eta << endl;

		// update wtx <- wtx + eta*qtx
		for (size_t c=0; c<numClasses; c++)
			for (size_t i=0; i<wtx[c].size(); i++)
				wtx[c][i]+=eta*qtx[c][i];

		// update wOld<- w ; w <- w + eta *q ; gOld<-g
		for (size_t c=0; c<numClasses; c++)
		{
			memcpy(&wOld[c][0],&w[c][0],sizeof(double)*w[c].size());
			memcpy(&gOld[c][0],&g[c][0],sizeof(double)*g[c].size());
			for (size_t i=0; i<numFeatures; i++)
				w[c][i]+= eta*q[c][i];
		}


		// initialize memory
		vector< vector< vector<double> > > memoryU(memorySize, vector< vector<double> >(numClasses));
		vector< vector< vector<double> > > memoryD(memorySize, vector< vector<double> >(numClasses));
		vector< double > memoryAlpha(memorySize);
		size_t nextMemPosition=0;
		size_t numMemPushes=0;
		
		// iterate until convergence
		size_t round=1;
		while (round<10000)
		{
			// compute errors and report round results
			{
				double trainingLogLikelihood=0.0, testLogLikelihood=0.0;
				const double trainingError = calcErrorRateWithLogLikelihood(trainingDataSet, trainingIdxs,
																		false, &trainingLogLikelihood);
				double testError=1.0;
				if (performTest)
					testError = calcErrorRateWithLogLikelihood(testDataSet, testIdxs, false, &testLogLikelihood);

				if (reportFrequency>0 && round % reportFrequency == 0)
				{
					cout << round << "\t" << scientific << setprecision(5) << trainingLogLikelihood << "\t" << fixed << setprecision(5) << trainingError;
					if (performTest)
						cout <<"\t" << scientific << testLogLikelihood << "\t" << fixed << setprecision(5)<< testError;
					cout << endl;
				}
				
				if (performTest)
				{
					if (testError<=bestTestError)
					{
						bestTestRound=round;
						bestTestError=testError;
						for (size_t c=0; c<numClasses; c++)
							memcpy(&bestW[c][0],&w[c][0],numFeatures*sizeof(double)); // copy weights
					}
				}
				
				if (trainingError<=bestTrainingError)
				{
					bestTrainingRound=round;
					bestTrainingError=trainingError;
					if (! performTest)
					{
						for (size_t c=0; c<numClasses; c++)
							memcpy(&bestW[c][0],&w[c][0],numFeatures*sizeof(double)); // copy weights
					}
				}		
			}

			// Train new round

			computeGradient(trainingDataSet, trainingIdxs, w, wtx, lambda, g);

			double alpha=0.0;
			double sigma=0.0;
			double utu=0.0;

			// write u=g'-g and d=w'-w onto memory, use them to compute alpha and sigma
			vector< vector<double> >& u = memoryU[nextMemPosition];
			vector< vector<double> >& d = memoryD[nextMemPosition];
			for (size_t c=0; c<numClasses; c++)
			{
				const size_t numFeatures = g[c].size();
				u[c].resize(numFeatures);
				d[c].resize(numFeatures);
				for (size_t i=0; i<numFeatures; i++)
				{
					const double gDiff = g[c][i]-gOld[c][i];
					const double wDiff = w[c][i]-wOld[c][i];
					u[c][i]=gDiff;
					d[c][i]=wDiff;
					alpha += gDiff*wDiff;
					utu += gDiff*gDiff;
				}
			}
			sigma = alpha / utu;
			memoryAlpha[nextMemPosition]=alpha;

			// update memory position
			nextMemPosition++;
			if (nextMemPosition == memorySize)
				nextMemPosition = 0;
			numMemPushes++;

			// q<-g
			for (size_t c=0; c<numClasses; c++)
				memcpy(&q[c][0],&g[c][0],g[c].size()*sizeof(double));
			
			// determine memory evaluation order 1..M (M is the newest)
			vector<size_t> memOrder;
			if (numMemPushes<=memorySize)
			{
				for (size_t i=0; i<numMemPushes; i++)
					memOrder.push_back(i);
			}
			else
			{
				for (size_t i=0; i<memorySize; i++)
					memOrder.push_back((i+nextMemPosition) % memorySize);
			}

			vector<double> beta(memOrder.size(),0.0);
			for (int i=memOrder.size()-1; i>=0; i--)
			{
				const size_t m = memOrder[static_cast<size_t>(i)];
				const double alpha = memoryAlpha[m];
				
				const vector< vector<double> >& dM = memoryD[m];
				double& betaM = beta[m];

				// compute beta[m] = (memory_d[m] dot g)/alpha[m]
				for (size_t c=0; c<dM.size(); c++)
					for (size_t i=0; i<dM[c].size(); i++)
						betaM += dM[c][i]*g[c][i];
				betaM/=alpha;
				
				// q <- q - beta[m]*memory_u[m]
				const vector< vector<double> >& uM = memoryU[m]; 
				for (size_t c=0; c<q.size(); c++)
					for (size_t i=0; i<q[c].size(); i++)
						q[c][i] -= betaM * uM[c][i];

			}

			// q <- sigma*q
			for (size_t c=0; c<q.size(); c++)
				for (size_t i=0; i<q[c].size(); i++)
					q[c][i]*=sigma;


			for (size_t i=0; i<memOrder.size(); i++)
			{
				const size_t m = memOrder[static_cast<size_t>(i)];
				const vector< vector<double> >& uM = memoryU[m];
				const vector< vector<double> >& dM = memoryD[m]; 
				const double betaM = beta[m];
				const double oneOverAlpha = 1.0 / memoryAlpha[m];
				double umq = computeDotProduct(uM,q);
				for (size_t c=0; c<numClasses; c++)
					for (size_t j=0; j<q[c].size(); j++)
					{
						const double dq = dM[c][j] * (betaM - umq*oneOverAlpha);
						umq += uM[c][j]*dq;
						q[c][j] += dq;
					}
		
			}

			// q<- -q
			for (size_t c=0; c<numClasses; c++)
				for (size_t i=0; i<q[c].size(); i++)
					q[c][i]=-q[c][i];
			
			// qtx = q*X
			for (size_t i=0; i<trainingIdxs.size(); i++)
			{
				const MlSample& sample = trainingDataSet->getSample(trainingIdxs[i]);
				for (size_t c=0; c<numClasses; c++)
					qtx[c][i]=computeDotProduct(q[c],sample.pairs);
			}
			

			bool needToRestart=false;
			eta = lineSearch(trainingDataSet, trainingIdxs, w, wtx, qtx, g, q, lambda);
			if (eta<= 0.0)
			{
				// restart ?
				needToRestart = true;
			}

			// update wOld<- w ; w <- w + eta *q ; gOld<- g
			for (size_t c=0; c<numClasses; c++)
			{
				memcpy(&wOld[c][0],&w[c][0],sizeof(double)*w[c].size());
				memcpy(&gOld[c][0],&g[c][0],sizeof(double)*g[c].size());
				for (size_t i=0; i<numFeatures; i++)
					w[c][i]+= eta*q[c][i];
			}

			for (size_t c=0; c<numClasses; c++)
				for (size_t i=0; i<numSamples; i++)
					wtx[c][i]+=eta*qtx[c][i];

			round++;
			totalRounds++;
			if (terminateTraining || needToRestart)
				break;
		}
		
		if (terminateTraining)
			break;
	}

	if (! params->getIndHadInternalError())
	{
		params->setIndNormalTermination(true);
	}
	else
		cout << "Warning: encountered mathemtical error while training!" << endl;

	weights_ = bestW;

	cout << "W=" << endl;
	printVector(weights_);
	cout << endl;

	cout << "Terminated after " << totalRounds << " rounds (" << megaRound << " restarts)" << endl;
	cout << "Best training error  " << fixed << setprecision(8) << bestTrainingError << " (round " << bestTrainingRound << ")" << endl;
	if (performTest)
	cout << "Best test error      "  << bestTestError     << " (round " << bestTestRound << ")" << endl;

	indWasInitialized_ = true;

	//this->calcErrorRateWithPerplexity(trainingDataSet, trainingIdxs, true, NULL);
}
Exemplo n.º 4
0
double lineSearch(const MlDataSet* ds,
				   const vector<size_t>& idxs,
				   const vector< vector<double> >& w,    // the point x
				   const vector< vector<double> >& wtx,
				   const vector< vector<double> >& qtx,
				   const vector< vector<double> >& g,    // the gradient
				   const vector< vector<double> >& q,    // q=p, direction along which we want to minimize f(x)
				   double lambda)
{
	const double TOLX = numeric_limits<double>::epsilon();
	const double ALF  = 1e-4;
	double wtw=0.0, qtw=0.0;
	double qtg = computeDotProduct(q,g);
	double qtq = computeDotProduct(q,q);
	double scale = 100.0/sqrt(qtq);
	if (scale>1.0)
		scale=1.0;

	qtq *= (scale * scale);
	double slope = scale * qtg;
	if (slope<0.0)
		return 0.0;

	if (lambda != 0.0)
	{
		wtw = computeDotProduct(w,w);
		qtw = computeDotProduct(q,w)*scale;
	}
	const double fOld  = computePosterior(ds, idxs, wtx, qtx, wtw, qtw, qtq, lambda, 0.0);
	
	double test=0.0;
	for (size_t c=0; c<q.size(); c++)
		for (size_t i=0; i<q[c].size(); i++)
		{
			double maxArg = fabs(w[c][i]);
			if (maxArg<1.0)
				maxArg=1.0;
			double temp=(scale*(fabs(q[c][i])/maxArg));
			if (temp>test)
				test=temp;
		}


	double alamMin = TOLX/test;
	double alam = 1.0;
	double alam2 = 0.0;
	double f2 = fOld;
	while (true)
	{
	//	double slopeAlam = slope * alam;
		if (alam<alamMin)
			return (alamMin*(scale<1.0 ? scale : 1.0));

		double f=computePosterior(ds, idxs, wtx, qtx, wtw, qtw, qtq, lambda, alam*scale);
		if (f>= fOld + ALF * alam * scale * slope)
		{
			return (alam*scale);
		}

		double tmpAlam;
		if (fabs(alam-1.0)<1e-20)
		{
			tmpAlam = -slope / (2.0 *(f - fOld - slope));
		}
		else
		{
			double r1 = f - fOld - alam*slope;
			double r2 = f2 - fOld - alam2*slope;
			double a = (r1/(alam*alam) - r2/(alam2*alam2))/(alam-alam2);
			double b = (-(alam2*r1)/(alam*alam)+(alam*r2)/(alam2*alam2))/(alam-alam2);
			if (fabs(a)<1e-22)
			{
				tmpAlam=-slope/(2.0*b);
			}
			else
			{
				double disc=b*b-3*a*slope;
				if (disc<0)
				{
					tmpAlam = 0.5 * alam;
				}
				else if (b<=0.0)
				{
					tmpAlam = (sqrt(disc)-b)/(3.0*a);
				}
				else
					tmpAlam = -slope / (b + sqrt(disc));
			}
			if (alam * 0.5 < tmpAlam)
				tmpAlam = alam *0.5;
		}
		
		alam2 = alam;
		f2 = f;
		alam = 0.1 * alam;
		if (alam<tmpAlam)
			alam = tmpAlam;
		
	}
	return (alam*scale);//
}