Пример #1
0
void PCM::propagateCavityGradients(const ScalarField& A_shape, ScalarField& A_nCavity, ScalarFieldTilde& A_rhoExplicitTilde, bool electricOnly) const
{   if(fsp.pcmVariant == PCM_SGA13)
    {   //Propagate gradient w.r.t expanded cavities to nCavity:
        ((PCM*)this)->A_nc = 0;
        const ScalarField* A_shapeEx[2] = { &A_shape, &Acavity_shapeVdw };
        for(int i=0; i<2; i++)
        {   //First compute derivative w.r.t expanded electron density:
            ScalarField A_nCavityEx;
            ShapeFunction::propagateGradient(nCavityEx[i], *(A_shapeEx[i]), A_nCavityEx, fsp.nc, fsp.sigma);
            ((PCM*)this)->A_nc += (-1./fsp.nc) * integral(A_nCavityEx*nCavityEx[i]);
            //then propagate to original electron density:
            ScalarField nCavityExUnused; //unused return value below
            ShapeFunction::expandDensity(wExpand[i], Rex[i], nCavity, nCavityExUnused, &A_nCavityEx, &A_nCavity);
        }
    }
    else if(fsp.pcmVariant == PCM_CANDLE)
    {   ScalarField A_nCavityEx;
        ScalarFieldTilde A_phiExt;
        double A_pCavity=0.;
        ShapeFunction::propagateGradient(nCavityEx[0], coulomb(Sf[0]*rhoExplicitTilde), I(wExpand[0]*J(A_shape)) + Acavity_shapeVdw,
                                         A_nCavityEx, A_phiExt, A_pCavity, fsp.nc, fsp.sigma, fsp.pCavity);
        A_nCavity += fsp.Ztot * I(Sf[0] * J(A_nCavityEx));
        if(!electricOnly) A_rhoExplicitTilde += coulomb(Sf[0]*A_phiExt);
        ((PCM*)this)->A_nc = (-1./fsp.nc) * integral(A_nCavityEx*nCavityEx[0]);
        ((PCM*)this)->A_eta_wDiel = integral(A_shape * I(wExpand[1]*J(shapeVdw)));
        ((PCM*)this)->A_pCavity = A_pCavity;
    }
    else if(isPCM_SCCS(fsp.pcmVariant))
    {   //Electrostatic and volumetric combinations via shape:
        ShapeFunctionSCCS::propagateGradient(nCavity, A_shape - fsp.cavityPressure, A_nCavity, fsp.rhoMin, fsp.rhoMax, epsBulk);
        //Add surface contributions:
        ScalarField shapePlus, shapeMinus;
        ShapeFunctionSCCS::compute(nCavity+(0.5*fsp.rhoDelta), shapePlus, fsp.rhoMin, fsp.rhoMax, epsBulk);
        ShapeFunctionSCCS::compute(nCavity-(0.5*fsp.rhoDelta), shapeMinus, fsp.rhoMin, fsp.rhoMax, epsBulk);
        VectorField Dn = gradient(nCavity);
        ScalarField DnLength = sqrt(lengthSquared(Dn));
        ScalarField A_shapeMinus = (fsp.cavityTension/fsp.rhoDelta) * DnLength;
        ScalarField A_DnLength = (fsp.cavityTension/fsp.rhoDelta) * (shapeMinus - shapePlus);
        A_nCavity -= divergence(Dn * (inv(DnLength) * A_DnLength));
        ShapeFunctionSCCS::propagateGradient(nCavity+(0.5*fsp.rhoDelta), -A_shapeMinus, A_nCavity, fsp.rhoMin, fsp.rhoMax, epsBulk);
        ShapeFunctionSCCS::propagateGradient(nCavity-(0.5*fsp.rhoDelta),  A_shapeMinus, A_nCavity, fsp.rhoMin, fsp.rhoMax, epsBulk);
    }
    else //All gradients are w.r.t the same shape function - propagate them to nCavity (which is defined as a density product for SaLSA)
    {   ShapeFunction::propagateGradient(nCavity, A_shape + Acavity_shape, A_nCavity, fsp.nc, fsp.sigma);
        ((PCM*)this)->A_nc = (-1./fsp.nc) * integral(A_nCavity*nCavity);
    }
}
Пример #2
0
double SaLSA::get_Adiel_and_grad_internal(ScalarFieldTilde& Adiel_rhoExplicitTilde, ScalarFieldTilde& Adiel_nCavityTilde, IonicGradient* extraForces, bool electricOnly) const
{
	EnergyComponents& Adiel = ((SaLSA*)this)->Adiel;
	const ScalarFieldTilde& phi = state; // that's what we solved for in minimize

	//First-order correct estimate of electrostatic energy:
	ScalarFieldTilde phiExt = coulomb(rhoExplicitTilde);
	Adiel["Electrostatic"] = -0.5*dot(phi, O(hessian(phi))) + dot(phi - 0.5*phiExt, O(rhoExplicitTilde));
	
	//Gradient w.r.t rhoExplicitTilde:
	Adiel_rhoExplicitTilde = phi - phiExt;

	//The "cavity" gradient is computed by chain rule via the gradient w.r.t to the shape function:
	const auto& solvent = fsp.solvents[0];
	ScalarField Adiel_shape; ScalarFieldArray Adiel_siteShape(solvent->molecule.sites.size());
	for(int r=rStart; r<rStop; r++)
	{	const MultipoleResponse& resp = *response[r];
		ScalarField& Adiel_s = resp.iSite<0 ? Adiel_shape : Adiel_siteShape[resp.iSite];
		if(resp.l>6) die("Angular momenta l > 6 not supported.\n");
		double prefac = 0.5 * 4*M_PI/(2*resp.l+1);
		ScalarFieldArray IlGradVphi = I(lGradient(resp.V * phi, resp.l));
		for(int lpm=0; lpm<(2*resp.l+1); lpm++)
			Adiel_s -= prefac * (IlGradVphi[lpm]*IlGradVphi[lpm]);
	}
	for(unsigned iSite=0; iSite<solvent->molecule.sites.size(); iSite++)
		if(Adiel_siteShape[iSite])
			Adiel_shape += I(Sf[iSite] * J(Adiel_siteShape[iSite]));
	nullToZero(Adiel_shape, gInfo); Adiel_shape->allReduce(MPIUtil::ReduceSum);
	
	//Propagate shape gradients to A_nCavity:
	ScalarField Adiel_nCavity;
	propagateCavityGradients(Adiel_shape, Adiel_nCavity, Adiel_rhoExplicitTilde, electricOnly);
	Adiel_nCavityTilde = nFluid * J(Adiel_nCavity);
	
	setExtraForces(extraForces, Adiel_nCavityTilde);
	return Adiel;
}
// elementary charge
inline quantity e() { return Rep( 1.602176462e-19L ) * coulomb(); }
double FluidMixture::operator()(const ScalarFieldArray& indep, ScalarFieldArray& Phi_indep, Outputs outputs) const
{	static StopWatch watch("FluidMixture::operator()"); watch.start();

	//logPrintf("indep.size: %d nIndep: %d\n",indep.size(),nIndep);
	myassert(indep.size()==get_nIndep());

	//---------- Compute site densities from the independent variables ---------
	ScalarFieldTildeArray Ntilde(nDensities); //site densities (in reciprocal space)
	std::vector< vector3<> > P0(component.size()); //polarization densities G=0
	for(unsigned ic=0; ic<component.size(); ic++)
	{	const FluidComponent& c = *component[ic];
		ScalarFieldArray N(c.molecule.sites.size());
		c.idealGas->getDensities(&indep[c.offsetIndep], &N[0], P0[ic]);
		for(unsigned i=0; i<c.molecule.sites.size(); i++)
		{	//Replace negative densities with 0:
			double Nmin, Nmax;
			callPref(eblas_capMinMax)(gInfo.nr, N[i]->dataPref(), Nmin, Nmax, 0.);
			//store site densities in fourier space
			Ntilde[c.offsetDensity+i] = J(N[i]);
			N[i] = 0; //Really skimp on memory!
		}
	}

	//----------- Handle density constraints ------------
	std::vector<double> Nscale(component.size(), 1.0); //density scale factor that satisfies the constraint
	std::vector<double> Ntot_c; //vector of total number of molecules for each component
	std::vector<double> Nscale_Qfixed(component.size(), 0.); //derivative of Nscale w.r.t the fixed charge
	std::vector<std::vector<double> > Nscale_N0(component.size(), std::vector<double>(component.size(),0.0)); //jacobian of Nscale w.r.t the uncorrected molecule counts
	std::vector<string> names; //list of molecule names
	
	//Find fixed N and charged species:
	double Qfixed = 0.0;
	if(rhoExternal) Qfixed += integral(rhoExternal);
	std::vector<std::pair<double,double> > N0Q;
	for(unsigned ic=0; ic<component.size(); ic++)
	{	const FluidComponent& c = *component[ic];
		double Qmolecule = c.molecule.getCharge();
		if(c.Nnorm>0 || Qmolecule)
		{	double N0 = integral(Ntilde[c.offsetDensity])/c.molecule.sites[0]->positions.size();
			if(c.Nnorm>0)
			{	Nscale[ic] = c.Nnorm/N0;
				Nscale_N0[ic][ic] = -c.Nnorm/pow(N0,2);
				Qfixed += Qmolecule*c.Nnorm;
			}
			else
			{	N0Q.push_back(std::make_pair(N0, Qmolecule));
				names.push_back(c.molecule.name);
			}
		}
	}
	//Find the betaV (see Qtot()) that makes the unit cell neutral
	if(N0Q.size()==0)
	{	if(fabs(Qfixed)>fabs(Qtol))
			die("Unit cell has a fixed net charge %le,"
				"and there are no free charged species to neutralize it.\n", Qfixed);
	}
	else
	{	double Qprime, Qcell, betaV=0.0;
		if(Qfixed+Qtot(-HUGE_VAL,Qprime,N0Q)<0)
			die("Unit cell will always have a net negative charge (no free positive charges).\n")
		if(Qfixed+Qtot(+HUGE_VAL,Qprime,N0Q)>0)
			die("Unit cell will always have a net positive charge (no free negative charges).\n")

		for(int iter=0; iter<10; iter++) //while(1)
		{	Qcell = Qfixed+Qtot(betaV,Qprime,N0Q,&names,&gInfo,verboseLog);
			if(verboseLog) logPrintf("betaV = %le, Qcell = %le, Qprime = %le\n", betaV, Qcell, Qprime);
			if(fabs(Qcell)<fabs(Qtol)) break;	
			if(std::isnan(Qcell)) 
			  { betaV=0.0; break;}
			betaV -= Qcell/Qprime; //Newton-Raphson update
			
		}
		for(unsigned ic=0; ic<component.size(); ic++)
		{	const FluidComponent& ci = *component[ic];
			double Qi = ci.molecule.getCharge();
			if(Qi && ci.Nnorm<=0)
			{	Nscale[ic] = exp(-Qi*betaV);
				Nscale_Qfixed[ic] = exp(-Qi*betaV) * Qi / Qprime;
				for(unsigned jc=0; jc<component.size(); jc++)
				{	const FluidComponent& cj = *component[jc];
					double Qj = cj.molecule.getCharge();
					if(Qj && cj.Nnorm<=0)
						Nscale_N0[ic][jc] += Qi*Qj*exp(-(Qi+Qj)*betaV)/Qprime;
				}
			}
		}
	}
	std::vector<double> Phi_Nscale(component.size(), 0.0); //accumulate explicit derivatives w.r.t Nscale here

	//Apply the scale factors to the site densities
	for(unsigned ic=0; ic<component.size(); ic++)
	{	const FluidComponent& c = *component[ic];
		for(unsigned i=0; i<c.molecule.sites.size(); i++)
			Ntilde[c.offsetDensity+i] *= Nscale[ic];
		P0[ic] *= Nscale[ic];
		Ntot_c.push_back(integral(Ntilde[c.offsetDensity]));
	}

	EnergyComponents Phi; //the grand free energy (with component information)
	ScalarFieldTildeArray Phi_Ntilde(nDensities); //gradients (functional derivative) w.r.t reciprocal space site densities
	nullToZero(Phi_Ntilde,gInfo);
	std::vector< vector3<> > Phi_P0(component.size()); //functional derivative w.r.t polarization density G=0
	VectorFieldTilde Phi_epsMF; //functional derivative w.r.t mean field electric field
	
	//G=0 fix for mismatch in fluid-fluid vs. fluid-electron charge kernels
	//We do this AFTER we have applied the appropriate scale factors
	if( N0Q.size() && (!useMFKernel)) //if we have charged species in our fluid, and are using mismatched charge kernels
	{
	  	for(unsigned ic=0; ic<component.size(); ic++)
		{	const FluidComponent& c = *component[ic];
		        for(unsigned i=0; i<c.molecule.sites.size(); i++)
			  { 
			    const Molecule::Site& s = *(c.molecule.sites[i]);
			    Phi["Gzero"] += Qfixed*(Ntot_c[ic]/gInfo.detR-c.idealGas->Nbulk)*s.positions.size()*s.deltaS;
			    Phi_Ntilde[c.offsetDensity+i]->data()[0] += (1.0/gInfo.dV) * (Qfixed*s.deltaS);
			  }
		}
	}

	//--------- Compute the (scaled) mean field coulomb interaction --------
	{	ScalarFieldTilde rho; //total charge density
		ScalarFieldTilde rhoMF; //effective charge density for mean-field term
		bool needRho = rhoExternal || outputs.Phi_rhoExternal;
		
		VectorFieldTilde epsMF = polarizable ? J(VectorField(&indep[nIndepIdgas])) : 0; //mean field electric field
		vector3<> P0tot;
		
		for(unsigned ic=0; ic<component.size(); ic++)
		{	const FluidComponent& c = *component[ic];
			for(unsigned i=0; i<c.molecule.sites.size(); i++)
			{	const Molecule::Site& s = *(c.molecule.sites[i]);
				if(s.chargeKernel)
				{	if(needRho) rho += s.chargeKernel * Ntilde[c.offsetDensity+i];
					rhoMF += s.chargeKernel(0) * (c.molecule.mfKernel * Ntilde[c.offsetDensity+i]);
				}
				//Polarization contributions:
				if(polarizable && s.polKernel)
				{	
					#define Polarization_Compute_Pi_Ni \
						VectorField Pi = (Cpol*s.alpha) * I(c.molecule.mfKernel*epsMF - (rhoExternal ? gradient(s.polKernel*coulomb(rhoExternal)) : 0)); \
						ScalarField Ni = I(Ntilde[c.offsetDensity+i]);
					Polarization_Compute_Pi_Ni
					
					ScalarField Phi_Ni = (0.5/(Cpol*s.alpha))*lengthSquared(Pi);
					Phi["Apol"] += gInfo.dV * dot(Ni, Phi_Ni);
					//Derivative contribution to site densities:
					Phi_Ntilde[c.offsetDensity+i] += Idag(Phi_Ni); Phi_Ni=0;
					//Update contributions to bound charge:
					VectorFieldTilde NPtilde = J(Ni * Pi); Pi=0; Ni=0;
					ScalarFieldTilde divNPbar;
					if(needRho) rho -= s.polKernel*divergence(NPtilde);
					VectorFieldTilde NPbarMF = c.molecule.mfKernel*NPtilde; NPtilde=0;
					rhoMF -= divergence(NPbarMF);
					Phi_epsMF += gInfo.nr * NPbarMF;
					P0tot += getGzero(NPbarMF);
				}
			}
			P0tot += P0[ic];
		}
		
		if(rhoMF)
		{	//External charge interaction:
			ScalarFieldTilde Phi_rho;
			ScalarFieldTilde Phi_rhoMF;
			if(needRho)
			{	if(rhoExternal)
				{	
				  ScalarFieldTilde OdExternal = O(coulomb(rhoExternal));
				  if (!useMFKernel)
				   {	
			  	       Phi["ExtCoulomb"] += dot(rho, OdExternal);
				       Phi_rho += OdExternal;
				   }
				   if (useMFKernel)
				   {
				       Phi["ExtCoulomb"] += dot(rhoMF, OdExternal);   
				       Phi_rhoMF += OdExternal; 
				   }
				}
				if(outputs.Phi_rhoExternal) 
				{
				  if (!useMFKernel) *outputs.Phi_rhoExternal = coulomb(rho);
				  if (useMFKernel) *outputs.Phi_rhoExternal = coulomb(rhoMF);
				}
			}
		
			//Mean field contributions:
			{	ScalarFieldTilde OdMF = O(coulomb(rhoMF)); //mean-field electrostatic potential
				Phi["Coulomb"] += 0.5*dot(rhoMF, OdMF);
				Phi_rhoMF += OdMF;
			}
			
			//Polarization density interactions:
			if(outputs.electricP) *outputs.electricP = P0tot * gInfo.detR;
			//--- corrections for net dipole in cell:
			vector3<> Phi_P0tot = (4*M_PI*gInfo.detR) * P0tot;
			Phi["PsqCell"] += 0.5 * dot(Phi_P0tot, P0tot); 
			//--- external electric field interactions:
			Phi["ExtCoulomb"] -= gInfo.detR * dot(Eexternal, P0tot);
			Phi_P0tot -= gInfo.detR * Eexternal;
			
			//Propagate gradients:
			for(unsigned ic=0; ic<component.size(); ic++)
			{	const FluidComponent& c = *component[ic];
				for(unsigned i=0; i<c.molecule.sites.size(); i++)
				{	const Molecule::Site& s = *(c.molecule.sites[i]);
					if(s.chargeKernel)
					{	if(Phi_rho) Phi_Ntilde[c.offsetDensity+i] += (1./gInfo.dV) * (s.chargeKernel * Phi_rho);
						Phi_Ntilde[c.offsetDensity+i] += (s.chargeKernel(0)/gInfo.dV) * (c.molecule.mfKernel * Phi_rhoMF);
					}
					//Polarization contributions:
					if(polarizable && s.polKernel)
					{	VectorFieldTilde Phi_NPtilde = gradient(c.molecule.mfKernel*Phi_rhoMF + (needRho ? s.polKernel*Phi_rho : 0));
						setGzero(Phi_NPtilde, getGzero(Phi_NPtilde) + Phi_P0tot);
						VectorField Phi_NP = Jdag(Phi_NPtilde); Phi_NPtilde=0;
						//propagate gradients from NP to N, epsMF and rhoExternal:
						Polarization_Compute_Pi_Ni
						#undef Polarization_Compute_Pi_Ni
						// --> via Ni
						ScalarField Phi_Ni; for(int k=0; k<3; k++) Phi_Ni += Phi_NP[k]*Pi[k];
						Phi_Ntilde[c.offsetDensity+i] += (1./gInfo.dV) * Idag(Phi_Ni); Phi_Ni=0;
						// --> via Pi
						VectorFieldTilde Phi_PiTilde = Idag(Phi_NP * Ni); Phi_NP=0;
						Phi_epsMF += (Cpol*s.alpha/gInfo.dV)*(c.molecule.mfKernel*Phi_PiTilde);
					}
				}
				Phi_P0[ic] += (1./gInfo.detR) * Phi_P0tot; //convert to functional derivative
			}
		}
	}
	
	//--------- Hard sphere mixture and bonding -------------
	{	//Compute the FMT weighted densities:
		ScalarFieldTilde n0tilde, n1tilde, n2tilde, n3tilde, n1vTilde, n2mTilde;
		std::vector<ScalarField> n0mol(component.size(), 0); //partial n0 for molecules that need bonding corrections
		std::vector<int> n0mult(component.size(), 0); //number of sites which contribute to n0 for each molecule
		std::vector<std::map<double,int> > bond(component.size()); //sets of bonds for each molecule
		bool bondsPresent = false; //whether bonds are present for any molecule
		for(unsigned ic=0; ic<component.size(); ic++)
		{	const FluidComponent& c = *component[ic];
			bond[ic] = c.molecule.getBonds();
			ScalarFieldTilde n0molTilde;
			for(unsigned i=0; i<c.molecule.sites.size(); i++)
			{	const Molecule::Site& s = *(c.molecule.sites[i]);
				if(s.Rhs)
				{	const ScalarFieldTilde& Nsite = Ntilde[c.offsetDensity+i];
					n0mult[ic] += s.positions.size();
					n0molTilde += s.w0  * Nsite;
					n1tilde    += s.w1  * Nsite;
					n2tilde    += s.w2  * Nsite;
					n3tilde    += s.w3  * Nsite;
					n1vTilde   += s.w1v * Nsite;
					n2mTilde   += s.w2m * Nsite;
				}
			}
			if(n0molTilde) n0tilde += n0molTilde;
			if(bond[ic].size())
			{	n0mol[ic] = I(n0molTilde);
				bondsPresent = true;
			}
		}
		if(n0tilde) //at least one sphere in the mixture
		{	ScalarField n0 = I(n0tilde); n0tilde=0;
			ScalarField n1 = I(n1tilde); n1tilde=0;
			ScalarField n2 = I(n2tilde); n2tilde=0;
			ScalarField Phi_n0, Phi_n1, Phi_n2; ScalarFieldTilde Phi_n3tilde, Phi_n1vTilde, Phi_n2mTilde;
			//Compute the sphere mixture free energy:
			Phi["MixedFMT"] += T * PhiFMT(n0, n1, n2, n3tilde, n1vTilde, n2mTilde,
				Phi_n0, Phi_n1, Phi_n2, Phi_n3tilde, Phi_n1vTilde, Phi_n2mTilde);
			//Bonding corrections if required
			if(bondsPresent)
			{	for(unsigned ic=0; ic<component.size(); ic++)
				{	const FluidComponent& c = *component[ic];
					ScalarField Phi_n0mol;
					for(const auto& b: bond[ic])
						Phi["Bonding"] += T * PhiBond(b.first, b.second*1./n0mult[ic],
							n0mol[ic], n2, n3tilde, Phi_n0mol, Phi_n2, Phi_n3tilde);
					if(Phi_n0mol)
					{	//Propagate gradient w.r.t n0mol[ic] to the site densities:
						ScalarFieldTilde Phi_n0molTilde = Idag(Phi_n0mol);
						for(unsigned i=0; i<c.molecule.sites.size(); i++)
						{	const Molecule::Site& s = *(c.molecule.sites[i]);
							if(s.Rhs)
								Phi_Ntilde[c.offsetDensity+i] += T * (s.w0  * Phi_n0molTilde);
						}
					}
				}
			}
			//Accumulate gradients w.r.t weighted densities to site densities:
			ScalarFieldTilde Phi_n0tilde = Idag(Phi_n0); Phi_n0=0;
			ScalarFieldTilde Phi_n1tilde = Idag(Phi_n1); Phi_n1=0;
			ScalarFieldTilde Phi_n2tilde = Idag(Phi_n2); Phi_n2=0;
			for(const FluidComponent* c: component)
			{	for(unsigned i=0; i<c->molecule.sites.size(); i++)
				{	const Molecule::Site& s = *(c->molecule.sites[i]);
					if(s.Rhs)
					{	ScalarFieldTilde& Phi_Nsite = Phi_Ntilde[c->offsetDensity+i];
						Phi_Nsite += T * (s.w0  * Phi_n0tilde);
						Phi_Nsite += T * (s.w1  * Phi_n1tilde);
						Phi_Nsite += T * (s.w2  * Phi_n2tilde);
						Phi_Nsite += T * (s.w3  * Phi_n3tilde);
						Phi_Nsite += T * (s.w1v * Phi_n1vTilde);
						Phi_Nsite += T * (s.w2m * Phi_n2mTilde);
					}
				}
			}
		}
	}

	//---------- Excess functionals --------------
	for(const FluidComponent* c: component) if(c->fex)
		Phi["Fex("+c->molecule.name+")"] += c->fex->compute(&Ntilde[c->offsetDensity], &Phi_Ntilde[c->offsetDensity]);

	//--------- Mixing functionals --------------
	for(const Fmix* fmix: fmixArr)
		Phi["Fmix("+fmix->getName()+")"] += fmix->compute(Ntilde, Phi_Ntilde);

	//--------- PhiNI ---------
	nullToZero(Phi_Ntilde, gInfo);
	if(outputs.N) outputs.N->resize(nDensities);
	//Put the site densities and gradients back in real space
	ScalarFieldArray N(nDensities);
	ScalarFieldArray Phi_N(nDensities);
	for(unsigned i=0; i<nDensities; i++)
	{	N[i] = I(Ntilde[i]); Ntilde[i]=0;
		Phi_N[i] = Jdag(Phi_Ntilde[i]); Phi_Ntilde[i] = 0;
		if(outputs.N) (*outputs.N)[i] = N[i]; //Link site-density to return pointer if necessary
	}
	//Estimate psiEff based on gradients, if requested
	if(outputs.psiEff)
	{	outputs.psiEff->resize(nDensities);
		for(const FluidComponent* c: component)
			for(unsigned i=0; i<c->molecule.sites.size(); i++)
			{	ScalarField& psiCur = outputs.psiEff->at(c->offsetDensity+i);
				psiCur = Phi_N[c->offsetDensity+i] + c->idealGas->V[i];
				if(i==0) psiCur -= c->idealGas->mu / c->molecule.sites[0]->positions.size();
				psiCur *= (-1./T);
			}
	}
	for(unsigned ic=0; ic<component.size(); ic++)
	{	const FluidComponent& c = *component[ic];
		Phi["PhiNI("+c.molecule.name+")"] +=
			c.idealGas->compute(&indep[c.offsetIndep], &N[c.offsetDensity], &Phi_N[c.offsetDensity], Nscale[ic], Phi_Nscale[ic]);

		//Fixed N correction to entropy:
		if(Nscale[ic]!=1.0)
		{	double deltaTs = T*log(Nscale[ic]) / c.molecule.sites[0]->positions.size();
			Phi_N[c.offsetDensity] += deltaTs;
			Phi["PhiNI("+c.molecule.name+")"] += integral(N[c.offsetDensity])*deltaTs;
		}
	}
	//Add in the implicit contributions to Phi_Nscale
	for(unsigned ic=0; ic<component.size(); ic++)
	{	const FluidComponent& ci = *component[ic];
		bool anyNonzero=false;
		for(unsigned jc=0; jc<component.size(); jc++)
			if(Nscale_N0[ic][jc])
				anyNonzero=true;
		if(anyNonzero)
		{	Phi_Nscale[ic] += gInfo.detR*dot(P0[ic], Phi_P0[ic])/ Nscale[ic];
			for(unsigned i=0; i<ci.molecule.sites.size(); i++)
				Phi_Nscale[ic] += gInfo.dV*dot(N[ci.offsetDensity+i], Phi_N[ci.offsetDensity+i])/ Nscale[ic];
		}
	}
	//Propagate gradients from Nscale to N:
	for(unsigned jc=0; jc<component.size(); jc++)
	{	const FluidComponent& cj = *component[jc];
		double Phi_Ncontrib = 0.0;
		for(unsigned ic=0; ic<component.size(); ic++)
			if(Nscale_N0[ic][jc])
				Phi_Ncontrib += Phi_Nscale[ic] * Nscale_N0[ic][jc];
		if(Phi_Ncontrib)
			Phi_N[cj.offsetDensity] += Phi_Ncontrib / (Nscale[jc] * cj.molecule.sites[0]->positions.size());
	}

	//Propagate gradients from Phi_N and Phi_P to Phi_indep
	Phi_indep.resize(get_nIndep());
	for(unsigned ic=0; ic<component.size(); ic++)
	{	const FluidComponent& c = *component[ic];
		c.idealGas->convertGradients(&indep[c.offsetIndep], &N[c.offsetDensity],
			&Phi_N[c.offsetDensity], Phi_P0[ic], &Phi_indep[c.offsetIndep], Nscale[ic]);
	}
	for(unsigned k=nIndepIdgas; k<get_nIndep(); k++) Phi_indep[k] = Jdag(Phi_epsMF[k-nIndepIdgas]);
	
	//Propagate gradients from Nscale to Qfixed / rhoExternal (Natural G=0 solution)
	if(outputs.Phi_rhoExternal)
	{	double Phi_Qfixed = 0.;
		for(unsigned ic=0; ic<component.size(); ic++)
		{
			Phi_Qfixed += Phi_Nscale[ic] * Nscale_Qfixed[ic];
			if (N0Q.size() && (!useMFKernel))
		        {
				const FluidComponent& c = *component[ic];
				for(unsigned i=0; i<c.molecule.sites.size(); i++)
				{ 
				 	//Correction to lambda from charge kernel mismatch
					  const Molecule::Site& s = *(c.molecule.sites[i]);
					  double lambda_s =  (Ntot_c[ic]/gInfo.detR-c.idealGas->Nbulk)*s.deltaS*s.positions.size();
					  if(verboseLog) logPrintf("Charge kernel mismatch correction for site %s of molecule %s: %lg\n",
					  s.name.c_str(),c.molecule.name.c_str(),lambda_s);
					  Phi_Qfixed += lambda_s;
				}
				 if(verboseLog) logPrintf("Total number of molecules of type %s: %lg\n",c.molecule.name.c_str(),Ntot_c[ic]);
		        }
		}
		nullToZero(*outputs.Phi_rhoExternal, gInfo);
		(*outputs.Phi_rhoExternal)->setGzero(Phi_Qfixed);
	}
	
	Phi["+pV"] += p * gInfo.detR; //background correction

	if(verboseLog) Phi.print(globalLog, true, "\t\t\t\t%15s = %25.16lf\n");
	if(outputs.Phi) *(outputs.Phi) = Phi;
	
	Phi_indep *= gInfo.dV; //convert functional derivative to partial derivative
	watch.stop();
	return Phi;
}
Пример #5
0
float coulomb_conf(int ires, int iconf, int jres, int jconf, PROT prot)
{
    float  e = 0.0;
    int    iatom, jatom;
    ATOM   *connect13[MAX_CONNECTED2],    *connect14[MAX_CONNECTED3];
    ATOM   *iatom_p, *jatom_p;
    int    connect13_res[MAX_CONNECTED2], connect14_res[MAX_CONNECTED3];
    int    iconnect, jconnect, kconnect, n_connect13, n_connect14;
    
    if (ires<0 || ires>prot.n_res-1) printf("   Error! coulomb_conf(): residue index out of range in protein\n"),exit(-1);
    if (jres<0 || jres>prot.n_res-1) printf("   Error! coulomb_conf(): residue index out of range in protein\n"),exit(-1);
    if (iconf<0 || iconf>prot.res[ires].n_conf-1) printf("   Error! coulomb_conf(): conformer index out of range in protein\n"),exit(-1);
    if (jconf<0 || jconf>prot.res[jres].n_conf-1) printf("   Error! coulomb_conf(): conformer index out of range in protein\n"),exit(-1);
    for (iatom=0; iatom<prot.res[ires].conf[iconf].n_atom; iatom++) {
        iatom_p = &prot.res[ires].conf[iconf].atom[iatom];
        if (!iatom_p->on) continue;
        
        n_connect13 = 0;
        n_connect14 = 0;
        memset(connect13,    0,MAX_CONNECTED2*sizeof(void *));
        memset(connect13_res,0,MAX_CONNECTED2*sizeof(int));
        memset(connect14,    0,MAX_CONNECTED3*sizeof(void *));
        memset(connect14_res,0,MAX_CONNECTED3*sizeof(int));
        
        for (iconnect = 0; iconnect < MAX_CONNECTED; iconnect++) {
            if (!iatom_p->connect12[iconnect]) break;
            n_connect13++;
            connect13[n_connect13-1] = iatom_p->connect12[iconnect];
            connect13_res[n_connect13-1] = iatom_p->connect12_res[iconnect];
            
            for (jconnect = 0; jconnect < MAX_CONNECTED; jconnect++) {
                if (!iatom_p->connect12[iconnect]->connect12[jconnect]) break;
                n_connect13++;
                connect13[n_connect13-1] = iatom_p->connect12[iconnect]->connect12[jconnect];
                connect13_res[n_connect13-1] = iatom_p->connect12[iconnect]->connect12_res[jconnect];

                for (kconnect = 0; kconnect < MAX_CONNECTED; kconnect++) {
                    if (!iatom_p->connect12[iconnect]->connect12[jconnect]->connect12[kconnect]) break;
                    n_connect14++;
                    connect14[n_connect14-1] = iatom_p->connect12[iconnect]->connect12[jconnect]->connect12[kconnect];
                    connect14_res[n_connect14-1] = iatom_p->connect12[iconnect]->connect12[jconnect]->connect12_res[kconnect];
                }
            }
        }
        
        for (jatom=0; jatom<prot.res[jres].conf[jconf].n_atom; jatom++) {
            jatom_p = &prot.res[jres].conf[jconf].atom[jatom];
            if (!jatom_p->on) continue;
            
            for (iconnect = 0; iconnect < n_connect13; iconnect++) {
                if (jres == connect13_res[iconnect] && !strcmp(jatom_p->name,connect13[iconnect]->name)) break;
            }
            if (iconnect < n_connect13) continue;

            for (iconnect = 0; iconnect < n_connect14; iconnect++) {
                if (jres == connect14_res[iconnect] && !strcmp(jatom_p->name,connect14[iconnect]->name)) break;
            }
            if (iconnect < n_connect13) {
                e += coulomb(*iatom_p, *jatom_p) * env.factor_14lj;
                continue;
            }
            
                e += coulomb(*iatom_p, *jatom_p);
                //printf("iatom %d,jatom %d\n",iatom,jatom);
        }
    }
    
    return e;
}
        void TransformerSPRINGEMBEDDER::SpringEmbed(Graph& G, vector<float> dimension_limits, IntermediateStepHandler* intermediatestephandler)
        {
            vector<VectorND> tractions;
            int dimensions = dimension_limits.size();

            unstressed_spring_length = dimension_limits[0];
            for(int i = 1; i < dimensions; i++)
                unstressed_spring_length *= dimension_limits[i];
            unstressed_spring_length /= G.NumberOfVertices();
            unstressed_spring_length = pow(unstressed_spring_length, 1.0f/dimensions) / (dimensions * 2);

            // init
            srand ( time(NULL) );
            for(VertexIterator a = G.BeginVertices(); a != G.EndVertices(); a++)
            {
                Coordinates coordinates;
                for(int i = 0; i < dimensions; i++)
                    coordinates.push_back(fmod(rand(), dimension_limits[i]));
                (*a)->SetCoordinates(coordinates);

                // should make sure that no two vertices have the same position
                tractions.push_back(VectorND(dimensions));
            }

            float max_movement;
            iteration = 0;
            do
            {
                // compute movement
                max_movement = 0;
                int i = 0;
                for(VertexIterator a = G.BeginVertices(); a != G.EndVertices(); a++, i++)
                {

                    VectorND traction = tractions[i] * friction;

                    // compute forces on a by the other vertices
                    VectorND aCoordinates((*a)->GetCoordinates(0));
                    for(VertexIterator b = G.BeginVertices(); b != G.EndVertices(); b++)
                    {
                        if(b==a)
                            continue;

                        // force on a by vertex b
                        VectorND bCoordinates((*b)->GetCoordinates(0));
                        traction += coulomb(aCoordinates, bCoordinates);
                        if((*a)->Adjacent(*b))
                            traction += hooke(aCoordinates, bCoordinates);
                    }

                    tractions[i] = traction;
                }


                // execute movement
                VertexIterator a = G.BeginVertices();
                for(int i = 0; a != G.EndVertices(); a++, i++)
                {
                    Coordinates OldCoordinates = (*a)->GetCoordinates(0);
                    Coordinates NewCoordinates(dimensions);
                    for(int j = 0; j < dimensions; j++)
                        NewCoordinates[j] = max(0.0f,min(dimension_limits[j],   OldCoordinates[j] + delta * tractions[i][j]   ));
                    (*a)->SetCoordinates( NewCoordinates );

                    // for the loop-condition
                    float current_movement = 0.0f;
                    for(int j = 0; j < dimensions; j++)
                        current_movement += (OldCoordinates[j] - NewCoordinates[j]) * (OldCoordinates[j] - NewCoordinates[j]);
                    if(sqrt(current_movement) > max_movement)
                        max_movement = sqrt(current_movement);
                }


                if(intermediatestephandler != NULL)
                    intermediatestephandler->Handle(&G);

                if(++iteration > nextincrease)
                {
                    movement_threshold += 0.01f;
                    nextincrease *= 4;
                }
                if(iteration > max_iterations)
                    break;

            }
            while(max_movement > movement_threshold*unstressed_spring_length);

        }
Пример #7
0
void PCM::updateCavity()
{
    //Cavities from expanded densities for SGA13 variant:
    if(fsp.pcmVariant == PCM_SGA13)
    {   ScalarField* shapeEx[2] = { &shape, &shapeVdw };
        for(int i=0; i<2; i++)
        {   ShapeFunction::expandDensity(wExpand[i], Rex[i], nCavity, nCavityEx[i]);
            ShapeFunction::compute(nCavityEx[i], *(shapeEx[i]), fsp.nc, fsp.sigma);
        }
    }
    else if(fsp.pcmVariant == PCM_CANDLE)
    {   nCavityEx[0] = fsp.Ztot * I(Sf[0] * J(nCavity));
        ShapeFunction::compute(nCavityEx[0], coulomb(Sf[0]*rhoExplicitTilde), shapeVdw,
                               fsp.nc, fsp.sigma, fsp.pCavity); //vdW cavity
        shape = I(wExpand[0] * J(shapeVdw)); //dielectric cavity
    }
    else if(isPCM_SCCS(fsp.pcmVariant))
        ShapeFunctionSCCS::compute(nCavity, shape, fsp.rhoMin, fsp.rhoMax, epsBulk);
    else //Compute directly from nCavity (which is a density product for SaLSA):
        ShapeFunction::compute(nCavity, shape, fsp.nc, fsp.sigma);

    //Compute and cache cavitation energy and gradients:
    const auto& solvent = fsp.solvents[0];
    switch(fsp.pcmVariant)
    {
    case PCM_SaLSA:
    case PCM_CANDLE:
    case PCM_SGA13:
    {   //Select relevant shape function:
        const ScalarFieldTilde sTilde = J(fsp.pcmVariant==PCM_SaLSA ? shape : shapeVdw);
        ScalarFieldTilde A_sTilde;
        //Cavitation:
        const double nlT = solvent->Nbulk * fsp.T;
        const double Gamma = log(nlT/solvent->Pvap) - 1.;
        const double Cp = 15. * (solvent->sigmaBulk/(2*solvent->Rvdw * nlT) - (1+Gamma)/6);
        const double coeff2 = 1. + Cp - 2.*Gamma;
        const double coeff3 = Gamma - 1. -2.*Cp;
        ScalarField sbar = I(wCavity*sTilde);
        Adiel["Cavitation"] = nlT * integral(sbar*(Gamma + sbar*(coeff2 + sbar*(coeff3 + sbar*Cp))));
        A_sTilde += wCavity*Idag(nlT * (Gamma + sbar*(2.*coeff2 + sbar*(3.*coeff3 + sbar*(4.*Cp)))));
        //Dispersion:
        ScalarFieldTildeArray Ntilde(Sf.size()), A_Ntilde(Sf.size()); //effective nuclear densities in spherical-averaged ansatz
        for(unsigned i=0; i<Sf.size(); i++)
            Ntilde[i] = solvent->Nbulk * (Sf[i] * sTilde);
        const double vdwScaleEff = (fsp.pcmVariant==PCM_CANDLE) ? fsp.sqrtC6eff : fsp.vdwScale;
        Adiel["Dispersion"] = e.vanDerWaals->energyAndGrad(atpos, Ntilde, atomicNumbers, vdwScaleEff, &A_Ntilde);
        A_vdwScale = Adiel["Dispersion"]/vdwScaleEff;
        for(unsigned i=0; i<Sf.size(); i++)
            if(A_Ntilde[i])
                A_sTilde += solvent->Nbulk * (Sf[i] * A_Ntilde[i]);
        //Propagate gradients to appropriate shape function:
        (fsp.pcmVariant==PCM_SaLSA ? Acavity_shape : Acavity_shapeVdw) = Jdag(A_sTilde);
        break;
    }
    case PCM_GLSSA13:
    {   VectorField Dshape = gradient(shape);
        ScalarField surfaceDensity = sqrt(lengthSquared(Dshape));
        ScalarField invSurfaceDensity = inv(surfaceDensity);
        A_tension = integral(surfaceDensity);
        Adiel["CavityTension"] = A_tension * fsp.cavityTension;
        Acavity_shape = (-fsp.cavityTension)*divergence(Dshape*invSurfaceDensity);
        break;
    }
    case PCM_LA12:
    case PCM_PRA05:
        break; //no contribution
case_PCM_SCCS_any:
        {   //Volume contribution:
            Adiel["CavityPressure"] = fsp.cavityPressure * (gInfo.detR - integral(shape));
            //Surface contribution:
            ScalarField shapePlus, shapeMinus;
            ShapeFunctionSCCS::compute(nCavity+(0.5*fsp.rhoDelta), shapePlus, fsp.rhoMin, fsp.rhoMax, epsBulk);
            ShapeFunctionSCCS::compute(nCavity-(0.5*fsp.rhoDelta), shapeMinus, fsp.rhoMin, fsp.rhoMax, epsBulk);
            ScalarField DnLength = sqrt(lengthSquared(gradient(nCavity)));
            Adiel["CavityTension"] = (fsp.cavityTension/fsp.rhoDelta) * integral(DnLength * (shapeMinus - shapePlus));
            break;
        }
    }
}