/*************************************************************************************
 *
 *  flattenedmodel.h:  suitable for deriving flattened models
 *
 *  Copyright (c) 2007 Jakub Schwarzmeier
 *  University of West Bohemia, Department of General Physics
 *
 *  This program is free software; you can redistribute it and/or modify
 *  it under the terms of the GNU General Public License as published by
 *  the Free Software Foundation; either version 2 of the License, or (at
 *  your option) any later version.
 *
 *  This program is distributed in the hope that it will be useful, but
 *  WITHOUT ANY WARRANTY; without even the implied warranty of
 *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
 *  General Public License for more details.
 *
 *  You should have received a copy of the GNU General Public License
 *  along with this program; if not, write to the Free Software
 *  Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
 *
 ************************************************************************/

#include "../includes.h"

#include "../../inc/ics/flattenedmodel.h"

CFlattenedModel::CFlattenedModel()
{
	m_name = "FlattenedModel";
	m_toomreQ = 0.0;		//  initially cold disk
}

void CFlattenedModel::Activate( double toomreQ )
{
	m_toomreQ = toomreQ;
}

double CFlattenedModel::PrepCumulativeMass( double r )
{
	return 2.0 * M_PI * r * rho( r );
}

//  axisymmetric force, d(Phi)/d(R)
double CFlattenedModel::dphi_dr(double r)
{
	return (m_G * CumulativeMass( r )) / (r * r);
//	return gsl_spline_eval_deriv (m_Phir, r, m_accPhir);
}

void CFlattenedModel::RealizeModel()
{
	//  test energies
    double Epot = 0;
    double Ekin = 0;

	//  prepare uniform (normal, Gaussian) random number generator
	gsl_rng * rnd = gsl_rng_alloc ( gsl_rng_mt19937 );

	std::string file = StartBinary();

    double mass = m_totalMass / m_count;

	//  Calculate position and velocity of each star
	//  in cylindrical coordinates.
    for ( unsigned int i = 0; i < m_count; i++ )
    {
		double mx, radius;
		do
		{
			mx = gsl_rng_uniform_pos( rnd ) * m_totalMass; //  CumulativeMass( TIMES_RADII * m_rMax );
			radius = rm( mx );
		}
		while (( radius > (m_rCut) ) || ( radius < (m_rMax/100.0) ) );

		double phi = 0 + gsl_rng_uniform_pos( rnd ) * ( 2 * M_PI + 0.0 );

		double px = radius * cos( phi );
		double py = radius * sin( phi );
		double pz = 0;	

        double vx, vy, vz;
        vx = vy = vz = 0;

        double pot = Potential( radius );
        Epot += pot;
		double psi = -pot;

		double v = sqrt( dphi_dr( radius ) * radius );
/*		if ( _isnan( v ) )
		{
			i--;
			continue;
		}
*/

		vx = - v * sin( phi );
		vy = v * cos( phi );
		vz = 0;

		//-- Velocity dispersions (fluctuations, deviations)
		if (m_toomreQ != 0.0)
		{
			double omega = v/radius;		//  angular velocity
			double kappa = sqrt(2.0)*omega;	//  epicyclic frequency

			//  Radial
			double sigmaRad = ( m_toomreQ * 3.36 * m_G * sigma(radius) ) / kappa;
//			if (m_G == 1.0)	sigmaRad *= 0.00000000088878;

			double sigmaRadX = sigmaRad * sin( phi );
			double sigmaRadY = sigmaRad * cos( phi );

			vx = CSphericalModel::RandomGaussian( vx, sigmaRadX );
			vy = CSphericalModel::RandomGaussian( vy, sigmaRadY );

			//  Tangential (azimuthal)
			double sigmaTan = ( m_toomreQ * 3.36 * m_G * sigma(radius) ) / (2.0 * omega );
//			if (m_G == 1.0)	sigmaTan *= 0.00000000088878;

			double sigmaTanX = sigmaTan * sin( phi );
			double sigmaTanY = sigmaTan * cos( phi );

			vx = CSphericalModel::RandomGaussian( vx, sigmaTanX );
			vy = CSphericalModel::RandomGaussian( vy, sigmaTanY );

			//  Vertical (z-component)
			vz = CSphericalModel::RandomGaussian( vz, 0.5 * sigmaTan );
		}

		Ekin += 0.5 * ( vx * vx + vy * vy + vz * vz );

		//--  Write ICs binary file
		BinaryWrite( mass );

		BinaryWrite( px );
		BinaryWrite( py );
		BinaryWrite( pz );

		BinaryWrite( vx );
		BinaryWrite( vy );
		BinaryWrite( vz );
    }

	StopBinary();

	xmlNode * pRootNode = StartXml();
	//  insert something into XML tree if needed
	StopXml(pRootNode);

	//  free random number generator object
	gsl_rng_free( rnd );

	std::cout << "Q = " << m_toomreQ << std::endl;

	PrintFinish( Ekin, Epot );

	std::cout << "Press enter key to finish...\a";
	getchar();
}

double CFlattenedModel::rho( double r )
{
	return sigma( r );
}