1 | //
|
---|
2 | // functional.cc
|
---|
3 | //
|
---|
4 | // Copyright (C) 1997 Limit Point Systems, Inc.
|
---|
5 | //
|
---|
6 | // Author: Curtis Janssen <cljanss@limitpt.com>
|
---|
7 | // Maintainer: LPS
|
---|
8 | //
|
---|
9 | // This file is part of the SC Toolkit.
|
---|
10 | //
|
---|
11 | // The SC Toolkit is free software; you can redistribute it and/or modify
|
---|
12 | // it under the terms of the GNU Library General Public License as published by
|
---|
13 | // the Free Software Foundation; either version 2, or (at your option)
|
---|
14 | // any later version.
|
---|
15 | //
|
---|
16 | // The SC Toolkit is distributed in the hope that it will be useful,
|
---|
17 | // but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
18 | // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
19 | // GNU Library General Public License for more details.
|
---|
20 | //
|
---|
21 | // You should have received a copy of the GNU Library General Public License
|
---|
22 | // along with the SC Toolkit; see the file COPYING.LIB. If not, write to
|
---|
23 | // the Free Software Foundation, 675 Mass Ave, Cambridge, MA 02139, USA.
|
---|
24 | //
|
---|
25 | // The U.S. Government is granted a limited license as per AL 91-7.
|
---|
26 | //
|
---|
27 |
|
---|
28 | #ifdef __GNUC__
|
---|
29 | #pragma implementation
|
---|
30 | #endif
|
---|
31 |
|
---|
32 | #include <util/misc/math.h>
|
---|
33 |
|
---|
34 | #include <util/misc/formio.h>
|
---|
35 | #include <util/state/stateio.h>
|
---|
36 | #include <chemistry/qc/dft/functional.h>
|
---|
37 |
|
---|
38 | using namespace std;
|
---|
39 | using namespace sc;
|
---|
40 |
|
---|
41 | #ifndef HAVE_ISNAN
|
---|
42 | #define isnan(x) ((x)!=(x))
|
---|
43 | #endif
|
---|
44 |
|
---|
45 | #define MIN_DENSITY 1.e-14
|
---|
46 | #define MIN_GAMMA 1.e-24
|
---|
47 | #define MIN_SQRTGAMMA 1.e-12
|
---|
48 | #define MAX_ZETA 1.-1.e-12
|
---|
49 | #define MIN_ZETA -(1.-1.e-12)
|
---|
50 |
|
---|
51 | ///////////////////////////////////////////////////////////////////////////
|
---|
52 | // utility functions
|
---|
53 |
|
---|
54 | inline static double
|
---|
55 | norm(double v[3])
|
---|
56 | {
|
---|
57 | double x,y,z;
|
---|
58 | return sqrt((x=v[0])*x + (y=v[1])*y + (z=v[2])*z);
|
---|
59 | }
|
---|
60 |
|
---|
61 | inline static double
|
---|
62 | dot(double v[3], double w[3])
|
---|
63 | {
|
---|
64 | return v[0]*w[0] + v[1]*w[1] + v[2]*w[2];
|
---|
65 | }
|
---|
66 |
|
---|
67 | ///////////////////////////////////////////////////////////////////////////
|
---|
68 | // PointInputData
|
---|
69 |
|
---|
70 | void
|
---|
71 | PointInputData::compute_derived(int spin_polarized,
|
---|
72 | int need_gradient,
|
---|
73 | int need_hessian)
|
---|
74 | {
|
---|
75 | a.rho_13 = pow(a.rho, 1.0/3.0);
|
---|
76 | if (need_gradient) {
|
---|
77 | a.gamma = dot(a.del_rho,a.del_rho);
|
---|
78 | }
|
---|
79 | if (need_hessian) {
|
---|
80 | a.lap_rho = a.hes_rho[XX] + a.hes_rho[YY] + a.hes_rho[ZZ];
|
---|
81 | }
|
---|
82 |
|
---|
83 |
|
---|
84 | if (spin_polarized) {
|
---|
85 | b.rho_13 = pow(b.rho, 1.0/3.0);
|
---|
86 | if (need_gradient) {
|
---|
87 | b.gamma = dot(b.del_rho,b.del_rho);
|
---|
88 | }
|
---|
89 | if (need_hessian) {
|
---|
90 | b.lap_rho = b.hes_rho[XX] + b.hes_rho[YY] + b.hes_rho[ZZ];
|
---|
91 | }
|
---|
92 | }
|
---|
93 | else {
|
---|
94 | b = a;
|
---|
95 | if (need_gradient) gamma_ab = a.gamma;
|
---|
96 | }
|
---|
97 |
|
---|
98 | if (spin_polarized && need_gradient) {
|
---|
99 | gamma_ab = a.del_rho[0]*b.del_rho[0]
|
---|
100 | + a.del_rho[1]*b.del_rho[1]
|
---|
101 | + a.del_rho[2]*b.del_rho[2];
|
---|
102 | }
|
---|
103 |
|
---|
104 | }
|
---|
105 |
|
---|
106 |
|
---|
107 | ///////////////////////////////////////////////////////////////////////////
|
---|
108 | // DenFunctional
|
---|
109 |
|
---|
110 | static ClassDesc DenFunctional_cd(
|
---|
111 | typeid(DenFunctional),"DenFunctional",1,"public SavableState",
|
---|
112 | 0, 0, 0);
|
---|
113 |
|
---|
114 | DenFunctional::DenFunctional(StateIn& s):
|
---|
115 | SavableState(s)
|
---|
116 | {
|
---|
117 | s.get(a0_);
|
---|
118 | s.get(spin_polarized_);
|
---|
119 | s.get(compute_potential_);
|
---|
120 | }
|
---|
121 |
|
---|
122 | DenFunctional::DenFunctional()
|
---|
123 | {
|
---|
124 | a0_ = 0;
|
---|
125 | spin_polarized_ = 0;
|
---|
126 | compute_potential_ = 0;
|
---|
127 | }
|
---|
128 |
|
---|
129 | DenFunctional::DenFunctional(const Ref<KeyVal>& keyval)
|
---|
130 | {
|
---|
131 | // a0 is usually zero, except for ACM functionals.
|
---|
132 | a0_ = keyval->doublevalue("a0");
|
---|
133 | spin_polarized_ = 0;
|
---|
134 | compute_potential_ = 0;
|
---|
135 | }
|
---|
136 |
|
---|
137 | DenFunctional::~DenFunctional()
|
---|
138 | {
|
---|
139 | }
|
---|
140 |
|
---|
141 | void
|
---|
142 | DenFunctional::save_data_state(StateOut& s)
|
---|
143 | {
|
---|
144 | s.put(a0_);
|
---|
145 | s.put(spin_polarized_);
|
---|
146 | s.put(compute_potential_);
|
---|
147 | }
|
---|
148 |
|
---|
149 | double
|
---|
150 | DenFunctional::a0() const
|
---|
151 | {
|
---|
152 | return a0_;
|
---|
153 | }
|
---|
154 |
|
---|
155 | int
|
---|
156 | DenFunctional::need_density_gradient()
|
---|
157 | {
|
---|
158 | return 0;
|
---|
159 | }
|
---|
160 |
|
---|
161 | int
|
---|
162 | DenFunctional::need_density_hessian()
|
---|
163 | {
|
---|
164 | return 0;
|
---|
165 | }
|
---|
166 |
|
---|
167 | void
|
---|
168 | DenFunctional::set_spin_polarized(int i)
|
---|
169 | {
|
---|
170 | spin_polarized_ = i;
|
---|
171 | }
|
---|
172 |
|
---|
173 | void
|
---|
174 | DenFunctional::set_compute_potential(int i)
|
---|
175 | {
|
---|
176 | compute_potential_ = i;
|
---|
177 | }
|
---|
178 |
|
---|
179 | void
|
---|
180 | DenFunctional::gradient(const PointInputData& id, PointOutputData& od,
|
---|
181 | double *grad_f, int acenter,
|
---|
182 | GaussianBasisSet *basis,
|
---|
183 | const double *dmat_a, const double *dmat_b,
|
---|
184 | int ncontrib, const int *contrib,
|
---|
185 | int ncontrib_bf, const int *contrib_bf,
|
---|
186 | const double *bs, const double *bsg,
|
---|
187 | const double *bsh)
|
---|
188 | {
|
---|
189 | int need_gamma_terms = need_density_gradient();
|
---|
190 | point(id, od);
|
---|
191 | memset(grad_f, 0, sizeof(double)*basis->ncenter()*3);
|
---|
192 | #if 0
|
---|
193 | ExEnv::outn() << scprintf("gradient: rho_a= %12.8f rho_b= %12.8f need_gamma = %d",
|
---|
194 | id.a.rho, id.b.rho, need_gamma_terms) << endl;
|
---|
195 | ExEnv::outn() << scprintf(" gamma_aa= %12.8f gamma_bb= %12.8f gamma_ab= % 12.8f",
|
---|
196 | id.a.gamma, id.b.gamma, id.gamma_ab) << endl;
|
---|
197 | ExEnv::outn() << scprintf(" df_drho_a= % 12.8f df_drho_b= % 12.8f",
|
---|
198 | od.df_drho_a, od.df_drho_b) << endl;
|
---|
199 | ExEnv::outn() << scprintf(" df_dg_aa= % 12.8f df_dg_bb= % 12.8f df_dg_ab= % 12.8f",
|
---|
200 | od.df_dgamma_aa,od.df_dgamma_bb,od.df_dgamma_ab) << endl;
|
---|
201 | #endif
|
---|
202 |
|
---|
203 | if (need_gamma_terms) {
|
---|
204 | double drhoa = od.df_drho_a;
|
---|
205 | double drhob = od.df_drho_b;
|
---|
206 | for (int nu=0; nu<ncontrib_bf; nu++) {
|
---|
207 | int nut = contrib_bf[nu];
|
---|
208 | int nuatom = basis->shell_to_center(basis->function_to_shell(nut));
|
---|
209 | double dfa_phi_nu = drhoa * bs[nu];
|
---|
210 | double dfb_phi_nu = drhob * bs[nu];
|
---|
211 | for (int mu=0; mu<ncontrib_bf; mu++) {
|
---|
212 | int mut = contrib_bf[mu];
|
---|
213 | int muatom
|
---|
214 | = basis->shell_to_center(basis->function_to_shell(mut));
|
---|
215 | if (muatom!=acenter) {
|
---|
216 | int nutmut
|
---|
217 | = (nut>mut?((nut*(nut+1))/2+mut):((mut*(mut+1))/2+nut));
|
---|
218 | double rho_a = dmat_a[nutmut];
|
---|
219 | double rho_b = dmat_b[nutmut];
|
---|
220 | int ixyz;
|
---|
221 | for (ixyz=0; ixyz<3; ixyz++) {
|
---|
222 | double contrib = -2.0*bsg[mu*3+ixyz]
|
---|
223 | * (rho_a*dfa_phi_nu + rho_b*dfb_phi_nu);
|
---|
224 | #define hoff(i,j) ((j)<(i)?((i)*((i)+1))/2+(j):((j)*((j)+1))/2+(i))
|
---|
225 | // gamma_aa contrib
|
---|
226 | if (need_gamma_terms) {
|
---|
227 | contrib
|
---|
228 | += 4.0 * od.df_dgamma_aa * rho_a
|
---|
229 | * ( - bsg[mu*3+ixyz]
|
---|
230 | * ( id.a.del_rho[0]*bsg[nu*3+0]
|
---|
231 | +id.a.del_rho[1]*bsg[nu*3+1]
|
---|
232 | +id.a.del_rho[2]*bsg[nu*3+2])
|
---|
233 | - bs[nu]
|
---|
234 | * ( id.a.del_rho[0]*bsh[mu*6+hoff(0,ixyz)]
|
---|
235 | +id.a.del_rho[1]*bsh[mu*6+hoff(1,ixyz)]
|
---|
236 | +id.a.del_rho[2]*bsh[mu*6+hoff(2,ixyz)]
|
---|
237 | )
|
---|
238 | );
|
---|
239 | }
|
---|
240 | // gamma_ab contrib
|
---|
241 | if (need_gamma_terms) {
|
---|
242 | contrib
|
---|
243 | += 2.0 * od.df_dgamma_ab * rho_a
|
---|
244 | * ( - bsg[mu*3+ixyz]
|
---|
245 | * ( id.b.del_rho[0]*bsg[nu*3+0]
|
---|
246 | +id.b.del_rho[1]*bsg[nu*3+1]
|
---|
247 | +id.b.del_rho[2]*bsg[nu*3+2])
|
---|
248 | - bs[nu]
|
---|
249 | * ( id.b.del_rho[0]*bsh[mu*6+hoff(0,ixyz)]
|
---|
250 | +id.b.del_rho[1]*bsh[mu*6+hoff(1,ixyz)]
|
---|
251 | +id.b.del_rho[2]*bsh[mu*6+hoff(2,ixyz)]
|
---|
252 | )
|
---|
253 | );
|
---|
254 | contrib
|
---|
255 | += 2.0 * od.df_dgamma_ab * rho_b
|
---|
256 | * ( - bsg[mu*3+ixyz]
|
---|
257 | * ( id.a.del_rho[0]*bsg[nu*3+0]
|
---|
258 | +id.a.del_rho[1]*bsg[nu*3+1]
|
---|
259 | +id.a.del_rho[2]*bsg[nu*3+2])
|
---|
260 | - bs[nu]
|
---|
261 | * ( id.a.del_rho[0]*bsh[mu*6+hoff(0,ixyz)]
|
---|
262 | +id.a.del_rho[1]*bsh[mu*6+hoff(1,ixyz)]
|
---|
263 | +id.a.del_rho[2]*bsh[mu*6+hoff(2,ixyz)]
|
---|
264 | )
|
---|
265 | );
|
---|
266 | }
|
---|
267 | // gamma_bb contrib
|
---|
268 | if (need_gamma_terms) {
|
---|
269 | contrib
|
---|
270 | += 4.0 * od.df_dgamma_bb * rho_b
|
---|
271 | * ( - bsg[mu*3+ixyz]
|
---|
272 | * ( id.b.del_rho[0]*bsg[nu*3+0]
|
---|
273 | +id.b.del_rho[1]*bsg[nu*3+1]
|
---|
274 | +id.b.del_rho[2]*bsg[nu*3+2])
|
---|
275 | - bs[nu]
|
---|
276 | * ( id.b.del_rho[0]*bsh[mu*6+hoff(0,ixyz)]
|
---|
277 | +id.b.del_rho[1]*bsh[mu*6+hoff(1,ixyz)]
|
---|
278 | +id.b.del_rho[2]*bsh[mu*6+hoff(2,ixyz)]
|
---|
279 | )
|
---|
280 | );
|
---|
281 | }
|
---|
282 | grad_f[3*muatom+ixyz] += contrib;
|
---|
283 | grad_f[3*acenter+ixyz] -= contrib;
|
---|
284 | }
|
---|
285 | }
|
---|
286 | }
|
---|
287 | }
|
---|
288 | }
|
---|
289 | else {
|
---|
290 | double drhoa = od.df_drho_a;
|
---|
291 | double drhob = od.df_drho_b;
|
---|
292 | for (int nu=0; nu<ncontrib_bf; nu++) {
|
---|
293 | int nut = contrib_bf[nu];
|
---|
294 | int nuatom = basis->shell_to_center(basis->function_to_shell(nut));
|
---|
295 | double dfa_phi_nu = drhoa * bs[nu];
|
---|
296 | double dfb_phi_nu = drhob * bs[nu];
|
---|
297 | for (int mu=0; mu<ncontrib_bf; mu++) {
|
---|
298 | int mut = contrib_bf[mu];
|
---|
299 | int muatom
|
---|
300 | = basis->shell_to_center(basis->function_to_shell(mut));
|
---|
301 | if (muatom!=acenter) {
|
---|
302 | int nutmut
|
---|
303 | = (nut>mut?((nut*(nut+1))/2+mut):((mut*(mut+1))/2+nut));
|
---|
304 | double rho_a = dmat_a[nutmut];
|
---|
305 | double rho_b = dmat_b[nutmut];
|
---|
306 | int ixyz;
|
---|
307 | for (ixyz=0; ixyz<3; ixyz++) {
|
---|
308 | // std::cout << "bsg[mu*3+ixyz] = " << bsg[mu*3+ixyz]
|
---|
309 | // << std::endl;
|
---|
310 | // std::cout << "rho_a = " << rho_a
|
---|
311 | // << std::endl;
|
---|
312 | // std::cout << "rho_b = " << rho_b
|
---|
313 | // << std::endl;
|
---|
314 | // std::cout << "dfa_phi_nu = " << dfa_phi_nu
|
---|
315 | // << std::endl;
|
---|
316 | // std::cout << "dfb_phi_nu = " << dfb_phi_nu
|
---|
317 | // << std::endl;
|
---|
318 | double contrib = -2.0*bsg[mu*3+ixyz]
|
---|
319 | * (rho_a*dfa_phi_nu + rho_b*dfb_phi_nu);
|
---|
320 | grad_f[3*muatom+ixyz] += contrib;
|
---|
321 | grad_f[3*acenter+ixyz] -= contrib;
|
---|
322 | }
|
---|
323 | }
|
---|
324 | }
|
---|
325 | }
|
---|
326 | }
|
---|
327 | }
|
---|
328 |
|
---|
329 | void
|
---|
330 | DenFunctional::do_fd_point(PointInputData&id,
|
---|
331 | double&in,double&out,
|
---|
332 | double lower_bound, double upper_bound)
|
---|
333 | {
|
---|
334 | double delta = 0.0000000001;
|
---|
335 | PointOutputData tod;
|
---|
336 | double insave = in;
|
---|
337 |
|
---|
338 | point(id,tod);
|
---|
339 | double outsave = tod.energy;
|
---|
340 |
|
---|
341 | int spin_polarized_save = spin_polarized_;
|
---|
342 | set_spin_polarized(1);
|
---|
343 |
|
---|
344 | if (insave-delta>=lower_bound && insave+delta<=upper_bound) {
|
---|
345 | in = insave+delta;
|
---|
346 | id.compute_derived(1, need_density_gradient(), false);
|
---|
347 | point(id,tod);
|
---|
348 | double plus = tod.energy;
|
---|
349 |
|
---|
350 | in = insave-delta;
|
---|
351 | id.compute_derived(1, need_density_gradient(), false);
|
---|
352 | point(id,tod);
|
---|
353 | double minu = tod.energy;
|
---|
354 | out = 0.5*(plus-minu)/delta;
|
---|
355 | }
|
---|
356 | else if (insave+2*delta<=upper_bound) {
|
---|
357 | in = insave+delta;
|
---|
358 | id.compute_derived(1, need_density_gradient(), false);
|
---|
359 | point(id,tod);
|
---|
360 | double plus = tod.energy;
|
---|
361 |
|
---|
362 | in = insave+2*delta;
|
---|
363 | id.compute_derived(1, need_density_gradient(), false);
|
---|
364 | point(id,tod);
|
---|
365 | double plus2 = tod.energy;
|
---|
366 | out = 0.5*(4.0*plus-plus2-3.0*outsave)/delta;
|
---|
367 | }
|
---|
368 | else if (insave-2*delta>=lower_bound) {
|
---|
369 | in = insave-delta;
|
---|
370 | id.compute_derived(1, need_density_gradient(), false);
|
---|
371 | point(id,tod);
|
---|
372 | double minu = tod.energy;
|
---|
373 |
|
---|
374 | in = insave-2*delta;
|
---|
375 | id.compute_derived(1, need_density_gradient(), false);
|
---|
376 | point(id,tod);
|
---|
377 | double minu2 = tod.energy;
|
---|
378 | out = -0.5*(4.0*minu-minu2-3.0*outsave)/delta;
|
---|
379 | }
|
---|
380 | else {
|
---|
381 | // the derivative is not well defined for this case
|
---|
382 | out = -135711.;
|
---|
383 | }
|
---|
384 | in = insave;
|
---|
385 | id.compute_derived(1, need_density_gradient(), false);
|
---|
386 |
|
---|
387 | set_spin_polarized(spin_polarized_save);
|
---|
388 | }
|
---|
389 |
|
---|
390 | void
|
---|
391 | DenFunctional::fd_point(const PointInputData&id, PointOutputData&od)
|
---|
392 | {
|
---|
393 | PointInputData tid(id);
|
---|
394 |
|
---|
395 | // fill in the energy at the initial density values
|
---|
396 | point(id,od);
|
---|
397 |
|
---|
398 | ExEnv::out0() << scprintf("ra=%7.5f rb=%7.5f gaa=%7.5f gbb=%7.5f gab= % 9.7f",
|
---|
399 | id.a.rho, id.b.rho, id.a.gamma, id.b.gamma, id.gamma_ab)
|
---|
400 | << endl;
|
---|
401 |
|
---|
402 | double ga = tid.a.gamma;
|
---|
403 | double gb = tid.b.gamma;
|
---|
404 | double gab = tid.gamma_ab;
|
---|
405 | double sga = sqrt(ga);
|
---|
406 | double sgb = sqrt(gb);
|
---|
407 |
|
---|
408 | double g_a_lbound = -2*gab - gb;
|
---|
409 | if (gb > 0 && gab*gab/gb > g_a_lbound) g_a_lbound = gab*gab/gb;
|
---|
410 | if (g_a_lbound < 0) g_a_lbound = 0.0;
|
---|
411 |
|
---|
412 | double g_b_lbound = -2*gab - ga;
|
---|
413 | if (ga > 0 && gab*gab/ga > g_b_lbound) g_b_lbound = gab*gab/ga;
|
---|
414 | if (g_b_lbound < 0) g_b_lbound = 0.0;
|
---|
415 |
|
---|
416 | double g_ab_lbound = -0.5*(ga+gb);
|
---|
417 | if (-sga*sgb > g_ab_lbound) g_ab_lbound = -sga*sgb;
|
---|
418 | // if (-sga*sgb < g_ab_lbound) g_ab_lbound = -sga*sgb;
|
---|
419 | double g_ab_ubound = sga*sgb;
|
---|
420 |
|
---|
421 | do_fd_point(tid, tid.a.rho, od.df_drho_a, 0.0, 10.0);
|
---|
422 | do_fd_point(tid, tid.b.rho, od.df_drho_b, 0.0, 10.0);
|
---|
423 | do_fd_point(tid, tid.a.gamma, od.df_dgamma_aa, g_a_lbound, 10.0);
|
---|
424 | do_fd_point(tid, tid.b.gamma, od.df_dgamma_bb, g_b_lbound, 10.0);
|
---|
425 | do_fd_point(tid, tid.gamma_ab, od.df_dgamma_ab, g_ab_lbound, g_ab_ubound);
|
---|
426 | }
|
---|
427 |
|
---|
428 | static int
|
---|
429 | check(const char *name, double fd, double an, const char *class_name)
|
---|
430 | {
|
---|
431 | // -135711. flags an undefined FD
|
---|
432 | if (fd == -135711.) return 0;
|
---|
433 |
|
---|
434 | double err = fabs(fd - an);
|
---|
435 | ExEnv::out0() << scprintf("%20s: fd = % 12.8f an = % 12.8f", name, fd, an)
|
---|
436 | << endl;
|
---|
437 | if ((fabs(an) > 0.03 && err/fabs(an) > 0.03)
|
---|
438 | || ((fabs(an) <= 0.03) && err > 0.03)
|
---|
439 | #ifdef HAVE_ISNAN
|
---|
440 | || isnan(an)
|
---|
441 | #endif
|
---|
442 | ) {
|
---|
443 | ExEnv::out0() << scprintf("Error: %12s: fd = % 12.8f an = % 12.8f (%s)",
|
---|
444 | name, fd, an, class_name)
|
---|
445 | << endl;
|
---|
446 | return 1;
|
---|
447 | }
|
---|
448 | return 0;
|
---|
449 | }
|
---|
450 |
|
---|
451 | int
|
---|
452 | DenFunctional::test(const PointInputData &id)
|
---|
453 | {
|
---|
454 | PointOutputData fd_od;
|
---|
455 | fd_point(id,fd_od);
|
---|
456 | PointOutputData an_od;
|
---|
457 | point(id,an_od);
|
---|
458 | int r = 0;
|
---|
459 | r+=check("df_drho_a", fd_od.df_drho_a, an_od.df_drho_a, class_name());
|
---|
460 | r+=check("df_drho_b", fd_od.df_drho_b, an_od.df_drho_b, class_name());
|
---|
461 | r+=check("df_dgamma_aa",fd_od.df_dgamma_aa,an_od.df_dgamma_aa,class_name());
|
---|
462 | r+=check("df_dgamma_ab",fd_od.df_dgamma_ab,an_od.df_dgamma_ab,class_name());
|
---|
463 | r+=check("df_dgamma_bb",fd_od.df_dgamma_bb,an_od.df_dgamma_bb,class_name());
|
---|
464 | return r;
|
---|
465 | }
|
---|
466 |
|
---|
467 | int
|
---|
468 | DenFunctional::test()
|
---|
469 | {
|
---|
470 | int i, j, k, l, m;
|
---|
471 | set_compute_potential(1);
|
---|
472 | set_spin_polarized(0);
|
---|
473 | SCVector3 r = 0.0;
|
---|
474 | PointInputData id(r);
|
---|
475 |
|
---|
476 | for (i=0; i<6; i++) id.a.hes_rho[i] = 0.0;
|
---|
477 | id.a.lap_rho = 0.0;
|
---|
478 |
|
---|
479 | // del rho should not be used by any of the functionals
|
---|
480 | for (i=0; i<3; i++) id.a.del_rho[i] = id.b.del_rho[i] = 0.0;
|
---|
481 |
|
---|
482 | double testrho[] = { 0.0, 0.001, 0.5, -1 };
|
---|
483 | double testgamma[] = { 0.0, 0.001, 0.5, -1 };
|
---|
484 | double testgammaab[] = { -0.5, 0.0, 0.5, -1 };
|
---|
485 |
|
---|
486 | int ret = 0;
|
---|
487 |
|
---|
488 | ExEnv::out0() << "Testing with rho_a == rho_b" << endl;
|
---|
489 | for (i=0; testrho[i] != -1.0; i++) {
|
---|
490 | if (testrho[i] == 0.0) continue;
|
---|
491 | id.a.rho=testrho[i];
|
---|
492 | for (j=0; testgamma[j] != -1.0; j++) {
|
---|
493 | if (testgamma[j] > testrho[i]) continue;
|
---|
494 | id.a.gamma = testgamma[j];
|
---|
495 | id.compute_derived(0, need_density_gradient(), false);
|
---|
496 | ret += test(id);
|
---|
497 | }
|
---|
498 | }
|
---|
499 |
|
---|
500 | set_spin_polarized(1);
|
---|
501 | ExEnv::out0() << "Testing with rho_a != rho_b" << endl;
|
---|
502 | for (i=0; testrho[i] != -1.0; i++) {
|
---|
503 | id.a.rho=testrho[i];
|
---|
504 | for (j=0; testrho[j] != -1.0; j++) {
|
---|
505 | id.b.rho=testrho[j];
|
---|
506 | if (testrho[i]+testrho[j] == 0.0) continue;
|
---|
507 | for (k=0; testgamma[k] != -1.0; k++) {
|
---|
508 | if (testgamma[k] > testrho[i]) continue;
|
---|
509 | id.a.gamma = testgamma[k];
|
---|
510 | double sqrt_gamma_a = sqrt(id.a.gamma);
|
---|
511 | for (l=0; testgamma[l] != -1.0; l++) {
|
---|
512 | if (testgamma[l] > testrho[j]) continue;
|
---|
513 | id.b.gamma = testgamma[l];
|
---|
514 | double sqrt_gamma_b = sqrt(id.b.gamma);
|
---|
515 | for (m=0; testgammaab[m] != -1.0; m++) {
|
---|
516 | // constrain gamma_ab to values allowed by the
|
---|
517 | // current gamma_a and gamma_b
|
---|
518 | id.gamma_ab = testgammaab[m];
|
---|
519 | if (id.gamma_ab > sqrt_gamma_a*sqrt_gamma_b) {
|
---|
520 | id.gamma_ab = sqrt_gamma_a*sqrt_gamma_b;
|
---|
521 | }
|
---|
522 | if (id.gamma_ab < -0.5*(id.a.gamma+id.b.gamma)) {
|
---|
523 | id.gamma_ab = -0.5*(id.a.gamma+id.b.gamma);
|
---|
524 | }
|
---|
525 | if (id.gamma_ab < -sqrt_gamma_a*sqrt_gamma_b) {
|
---|
526 | id.gamma_ab = -sqrt_gamma_a*sqrt_gamma_b;
|
---|
527 | }
|
---|
528 | id.compute_derived(1, need_density_gradient(), false);
|
---|
529 | ret += test(id);
|
---|
530 | }
|
---|
531 | }
|
---|
532 | }
|
---|
533 | }
|
---|
534 | }
|
---|
535 | return ret;
|
---|
536 | }
|
---|
537 |
|
---|
538 | /////////////////////////////////////////////////////////////////////////////
|
---|
539 | // NElFunctional
|
---|
540 |
|
---|
541 | static ClassDesc NElFunctional_cd(
|
---|
542 | typeid(NElFunctional),"NElFunctional",1,"public DenFunctional",
|
---|
543 | 0, create<NElFunctional>, create<NElFunctional>);
|
---|
544 |
|
---|
545 | NElFunctional::NElFunctional(StateIn& s):
|
---|
546 | SavableState(s),
|
---|
547 | DenFunctional(s)
|
---|
548 | {
|
---|
549 | }
|
---|
550 |
|
---|
551 | NElFunctional::NElFunctional(const Ref<KeyVal>& keyval):
|
---|
552 | DenFunctional(keyval)
|
---|
553 | {
|
---|
554 | }
|
---|
555 |
|
---|
556 | NElFunctional::~NElFunctional()
|
---|
557 | {
|
---|
558 | }
|
---|
559 |
|
---|
560 | void
|
---|
561 | NElFunctional::save_data_state(StateOut& s)
|
---|
562 | {
|
---|
563 | DenFunctional::save_data_state(s);
|
---|
564 | }
|
---|
565 |
|
---|
566 | void
|
---|
567 | NElFunctional::point(const PointInputData &id,
|
---|
568 | PointOutputData &od)
|
---|
569 | {
|
---|
570 | od.zero();
|
---|
571 | od.energy = id.a.rho + id.b.rho;
|
---|
572 | }
|
---|
573 |
|
---|
574 | /////////////////////////////////////////////////////////////////////////////
|
---|
575 | // SumDenFunctional
|
---|
576 |
|
---|
577 | static ClassDesc SumDenFunctional_cd(
|
---|
578 | typeid(SumDenFunctional),"SumDenFunctional",1,"public DenFunctional",
|
---|
579 | 0, create<SumDenFunctional>, create<SumDenFunctional>);
|
---|
580 |
|
---|
581 | SumDenFunctional::SumDenFunctional(StateIn& s):
|
---|
582 | SavableState(s),
|
---|
583 | DenFunctional(s),
|
---|
584 | n_(0),
|
---|
585 | funcs_(0),
|
---|
586 | coefs_(0)
|
---|
587 | {
|
---|
588 | s.get(n_);
|
---|
589 | if (n_) {
|
---|
590 | s.get(coefs_);
|
---|
591 | funcs_ = new Ref<DenFunctional>[n_];
|
---|
592 | for (int i=0; i < n_; i++)
|
---|
593 | funcs_[i] << SavableState::restore_state(s);
|
---|
594 | }
|
---|
595 | }
|
---|
596 |
|
---|
597 | SumDenFunctional::SumDenFunctional() :
|
---|
598 | n_(0),
|
---|
599 | funcs_(0),
|
---|
600 | coefs_(0)
|
---|
601 | {
|
---|
602 | }
|
---|
603 |
|
---|
604 | SumDenFunctional::SumDenFunctional(const Ref<KeyVal>& keyval):
|
---|
605 | DenFunctional(keyval),
|
---|
606 | n_(0),
|
---|
607 | funcs_(0),
|
---|
608 | coefs_(0)
|
---|
609 | {
|
---|
610 | int ncoef = keyval->count("coefs");
|
---|
611 | int nfunc = keyval->count("funcs");
|
---|
612 | if (ncoef != nfunc && ncoef != 0) {
|
---|
613 | ExEnv::out0() << "SumDenFunctional: number of coefs and funcs differ" << endl;
|
---|
614 | abort();
|
---|
615 | }
|
---|
616 |
|
---|
617 | n_ = nfunc;
|
---|
618 | coefs_ = new double[n_];
|
---|
619 | funcs_ = new Ref<DenFunctional>[n_];
|
---|
620 | for (int i=0; i < n_; i++) {
|
---|
621 | if (ncoef)
|
---|
622 | coefs_[i] = keyval->doublevalue("coefs", i);
|
---|
623 | else
|
---|
624 | coefs_[i] = 1.0;
|
---|
625 | funcs_[i] << keyval->describedclassvalue("funcs", i);
|
---|
626 | }
|
---|
627 | }
|
---|
628 |
|
---|
629 | SumDenFunctional::~SumDenFunctional()
|
---|
630 | {
|
---|
631 | if (n_) {
|
---|
632 | for (int i=0; i < n_; i++) funcs_[i] = 0; // just in case
|
---|
633 | delete[] funcs_;
|
---|
634 | delete[] coefs_;
|
---|
635 | }
|
---|
636 | n_=0;
|
---|
637 | funcs_=0;
|
---|
638 | coefs_=0;
|
---|
639 | }
|
---|
640 |
|
---|
641 | void
|
---|
642 | SumDenFunctional::save_data_state(StateOut& s)
|
---|
643 | {
|
---|
644 | DenFunctional::save_data_state(s);
|
---|
645 | s.put(n_);
|
---|
646 | if (n_) {
|
---|
647 | s.put(coefs_, n_);
|
---|
648 | for (int i=0; i < n_; i++)
|
---|
649 | SavableState::save_state(funcs_[i].pointer(),s);
|
---|
650 | }
|
---|
651 | }
|
---|
652 |
|
---|
653 | double
|
---|
654 | SumDenFunctional::a0() const
|
---|
655 | {
|
---|
656 | double eff_a0 = a0_;
|
---|
657 | for (int i=0; i < n_; i++) {
|
---|
658 | eff_a0 += coefs_[i] * funcs_[i]->a0();
|
---|
659 | }
|
---|
660 | return eff_a0;
|
---|
661 | }
|
---|
662 |
|
---|
663 | int
|
---|
664 | SumDenFunctional::need_density_gradient()
|
---|
665 | {
|
---|
666 | for (int i=0; i < n_; i++)
|
---|
667 | if (funcs_[i]->need_density_gradient())
|
---|
668 | return 1;
|
---|
669 |
|
---|
670 | return 0;
|
---|
671 | }
|
---|
672 |
|
---|
673 | void
|
---|
674 | SumDenFunctional::set_spin_polarized(int p)
|
---|
675 | {
|
---|
676 | spin_polarized_ = p;
|
---|
677 | for (int i=0; i < n_; i++)
|
---|
678 | funcs_[i]->set_spin_polarized(p);
|
---|
679 | }
|
---|
680 |
|
---|
681 | void
|
---|
682 | SumDenFunctional::set_compute_potential(int val)
|
---|
683 | {
|
---|
684 | compute_potential_ = val;
|
---|
685 | for (int i=0; i < n_; i++)
|
---|
686 | funcs_[i]->set_compute_potential(val);
|
---|
687 | }
|
---|
688 |
|
---|
689 | void
|
---|
690 | SumDenFunctional::point(const PointInputData &id,
|
---|
691 | PointOutputData &od)
|
---|
692 | {
|
---|
693 | od.zero();
|
---|
694 | PointOutputData tmpod;
|
---|
695 | for (int i=0; i < n_; i++) {
|
---|
696 | funcs_[i]->point(id, tmpod);
|
---|
697 |
|
---|
698 | od.energy += coefs_[i] * tmpod.energy;
|
---|
699 | if (compute_potential_) {
|
---|
700 | od.df_drho_a += coefs_[i] * tmpod.df_drho_a;
|
---|
701 | od.df_drho_b += coefs_[i] * tmpod.df_drho_b;
|
---|
702 | od.df_dgamma_aa += coefs_[i] * tmpod.df_dgamma_aa;
|
---|
703 | od.df_dgamma_ab += coefs_[i] * tmpod.df_dgamma_ab;
|
---|
704 | od.df_dgamma_bb += coefs_[i] * tmpod.df_dgamma_bb;
|
---|
705 | }
|
---|
706 | }
|
---|
707 | }
|
---|
708 |
|
---|
709 | void
|
---|
710 | SumDenFunctional::print(ostream& o) const
|
---|
711 | {
|
---|
712 | o
|
---|
713 | << indent << "Sum of Functionals:" << endl;
|
---|
714 | o << incindent;
|
---|
715 | o << indent << scprintf("%+18.16f Hartree-Fock Exchange",a0_) << endl;
|
---|
716 | for (int i=0; i<n_; i++) {
|
---|
717 | o << indent << scprintf("%+18.16f",coefs_[i]) << endl;
|
---|
718 | o << incindent;
|
---|
719 | funcs_[i]->print(o);
|
---|
720 | o << decindent;
|
---|
721 | }
|
---|
722 | o << decindent;
|
---|
723 | }
|
---|
724 |
|
---|
725 | /////////////////////////////////////////////////////////////////////////////
|
---|
726 | // StdDenFunctional
|
---|
727 |
|
---|
728 | static ClassDesc StdDenFunctional_cd(
|
---|
729 | typeid(StdDenFunctional),"StdDenFunctional",1,"public SumDenFunctional",
|
---|
730 | 0, create<StdDenFunctional>, create<StdDenFunctional>);
|
---|
731 |
|
---|
732 | StdDenFunctional::StdDenFunctional(StateIn& s):
|
---|
733 | SavableState(s),
|
---|
734 | SumDenFunctional(s),
|
---|
735 | name_(0)
|
---|
736 | {
|
---|
737 | s.getstring(name_);
|
---|
738 | }
|
---|
739 |
|
---|
740 | StdDenFunctional::StdDenFunctional():
|
---|
741 | name_(0)
|
---|
742 | {
|
---|
743 | }
|
---|
744 |
|
---|
745 | void
|
---|
746 | StdDenFunctional::init_arrays(int n)
|
---|
747 | {
|
---|
748 | n_ = n;
|
---|
749 | funcs_ = new Ref<DenFunctional>[n_];
|
---|
750 | coefs_ = new double[n_];
|
---|
751 | for (int i=0; i<n_; i++) coefs_[i] = 1.0;
|
---|
752 | }
|
---|
753 |
|
---|
754 | StdDenFunctional::StdDenFunctional(const Ref<KeyVal>& keyval)
|
---|
755 | {
|
---|
756 | name_ = keyval->pcharvalue("name");
|
---|
757 | if (name_) {
|
---|
758 | if (!strcmp(name_,"HFK")) {
|
---|
759 | n_ = 0;
|
---|
760 | a0_ = 1.0;
|
---|
761 | }
|
---|
762 | else if (!strcmp(name_,"XALPHA")) {
|
---|
763 | init_arrays(1);
|
---|
764 | funcs_[0] = new XalphaFunctional;
|
---|
765 | }
|
---|
766 | else if (!strcmp(name_,"HFS")) {
|
---|
767 | init_arrays(1);
|
---|
768 | funcs_[0] = new SlaterXFunctional;
|
---|
769 | }
|
---|
770 | else if (!strcmp(name_,"HFB")) {
|
---|
771 | init_arrays(2);
|
---|
772 | funcs_[0] = new SlaterXFunctional;
|
---|
773 | funcs_[1] = new Becke88XFunctional;
|
---|
774 | }
|
---|
775 | else if (!strcmp(name_,"HFG96")) {
|
---|
776 | init_arrays(1);
|
---|
777 | funcs_[0] = new G96XFunctional;
|
---|
778 | }
|
---|
779 | else if (!strcmp(name_,"G96LYP")) {
|
---|
780 | init_arrays(2);
|
---|
781 | funcs_[0] = new G96XFunctional;
|
---|
782 | funcs_[1] = new LYPCFunctional;
|
---|
783 | }
|
---|
784 | else if (!strcmp(name_,"BLYP")) {
|
---|
785 | init_arrays(3);
|
---|
786 | funcs_[0] = new SlaterXFunctional;
|
---|
787 | funcs_[1] = new Becke88XFunctional;
|
---|
788 | funcs_[2] = new LYPCFunctional;
|
---|
789 | }
|
---|
790 | else if (!strcmp(name_,"SVWN1")) {
|
---|
791 | init_arrays(2);
|
---|
792 | funcs_[0] = new SlaterXFunctional;
|
---|
793 | funcs_[1] = new VWN1LCFunctional;
|
---|
794 | }
|
---|
795 | else if (!strcmp(name_,"SVWN1RPA")) {
|
---|
796 | init_arrays(2);
|
---|
797 | funcs_[0] = new SlaterXFunctional;
|
---|
798 | funcs_[1] = new VWN1LCFunctional(1);
|
---|
799 | }
|
---|
800 | else if (!strcmp(name_,"SVWN2")) {
|
---|
801 | init_arrays(2);
|
---|
802 | funcs_[0] = new SlaterXFunctional;
|
---|
803 | funcs_[1] = new VWN2LCFunctional;
|
---|
804 | }
|
---|
805 | else if (!strcmp(name_,"SVWN3")) {
|
---|
806 | init_arrays(2);
|
---|
807 | funcs_[0] = new SlaterXFunctional;
|
---|
808 | funcs_[1] = new VWN3LCFunctional;
|
---|
809 | }
|
---|
810 | else if (!strcmp(name_,"SVWN4")) {
|
---|
811 | init_arrays(2);
|
---|
812 | funcs_[0] = new SlaterXFunctional;
|
---|
813 | funcs_[1] = new VWN4LCFunctional;
|
---|
814 | }
|
---|
815 | else if (!strcmp(name_,"SVWN5")) {
|
---|
816 | init_arrays(2);
|
---|
817 | funcs_[0] = new SlaterXFunctional;
|
---|
818 | funcs_[1] = new VWN5LCFunctional;
|
---|
819 | }
|
---|
820 | else if (!strcmp(name_,"SPZ81")) {
|
---|
821 | init_arrays(2);
|
---|
822 | funcs_[0] = new SlaterXFunctional;
|
---|
823 | funcs_[1] = new PZ81LCFunctional;
|
---|
824 | }
|
---|
825 | else if (!strcmp(name_,"SPW92")) {
|
---|
826 | init_arrays(2);
|
---|
827 | funcs_[0] = new SlaterXFunctional;
|
---|
828 | funcs_[1] = new PW92LCFunctional;
|
---|
829 | }
|
---|
830 | else if (!strcmp(name_,"BPW91")) {
|
---|
831 | init_arrays(3);
|
---|
832 | funcs_[0] = new SlaterXFunctional;
|
---|
833 | funcs_[1] = new Becke88XFunctional;
|
---|
834 | funcs_[2] = new PW91CFunctional;
|
---|
835 | }
|
---|
836 | else if (!strcmp(name_,"BP86")) {
|
---|
837 | init_arrays(4);
|
---|
838 | funcs_[0] = new SlaterXFunctional;
|
---|
839 | funcs_[1] = new Becke88XFunctional;
|
---|
840 | funcs_[2] = new P86CFunctional;
|
---|
841 | funcs_[3] = new PZ81LCFunctional;
|
---|
842 | }
|
---|
843 | else if (!strcmp(name_,"B3LYP")) {
|
---|
844 | init_arrays(4);
|
---|
845 | a0_ = 0.2;
|
---|
846 | coefs_[0] = 0.8;
|
---|
847 | coefs_[1] = 0.72;
|
---|
848 | coefs_[2] = 0.19;
|
---|
849 | coefs_[3] = 0.81;
|
---|
850 | funcs_[0] = new SlaterXFunctional;
|
---|
851 | funcs_[1] = new Becke88XFunctional;
|
---|
852 | funcs_[2] = new VWN1LCFunctional(1);
|
---|
853 | funcs_[3] = new LYPCFunctional;
|
---|
854 | }
|
---|
855 | else if (!strcmp(name_,"KMLYP")) {
|
---|
856 | init_arrays(3);
|
---|
857 | a0_ = 0.557;
|
---|
858 | coefs_[0] = 0.443;
|
---|
859 | coefs_[1] = 0.552;
|
---|
860 | coefs_[2] = 0.448;
|
---|
861 | funcs_[0] = new SlaterXFunctional;
|
---|
862 | funcs_[1] = new VWN1LCFunctional(1);
|
---|
863 | funcs_[2] = new LYPCFunctional;
|
---|
864 | }
|
---|
865 | else if (!strcmp(name_,"B3PW91")) {
|
---|
866 | init_arrays(4);
|
---|
867 | a0_ = 0.2;
|
---|
868 | coefs_[0] = 0.8;
|
---|
869 | coefs_[1] = 0.72;
|
---|
870 | coefs_[2] = 0.81;
|
---|
871 | coefs_[3] = 0.19;
|
---|
872 | funcs_[0] = new SlaterXFunctional;
|
---|
873 | funcs_[1] = new Becke88XFunctional;
|
---|
874 | funcs_[2] = new PW91CFunctional;
|
---|
875 | funcs_[3] = new PW92LCFunctional;
|
---|
876 | }
|
---|
877 | else if (!strcmp(name_,"B3P86")) {
|
---|
878 | init_arrays(4);
|
---|
879 | a0_ = 0.2;
|
---|
880 | coefs_[0] = 0.8;
|
---|
881 | coefs_[1] = 0.72;
|
---|
882 | coefs_[2] = 0.81;
|
---|
883 | coefs_[3] = 1.0;
|
---|
884 | funcs_[0] = new SlaterXFunctional;
|
---|
885 | funcs_[1] = new Becke88XFunctional;
|
---|
886 | funcs_[2] = new P86CFunctional;
|
---|
887 | funcs_[3] = new VWN1LCFunctional(1);
|
---|
888 | }
|
---|
889 | else if (!strcmp(name_,"PBE")) {
|
---|
890 | init_arrays(2);
|
---|
891 | funcs_[0] = new PBEXFunctional;
|
---|
892 | funcs_[1] = new PBECFunctional;
|
---|
893 | }
|
---|
894 | else if (!strcmp(name_,"PW91")) {
|
---|
895 | init_arrays(2);
|
---|
896 | funcs_[0] = new PW91XFunctional;
|
---|
897 | funcs_[1] = new PW91CFunctional;
|
---|
898 | }
|
---|
899 | else if (!strcmp(name_,"mPW(PW91)PW91")) {
|
---|
900 | init_arrays(2);
|
---|
901 | funcs_[0] = new mPW91XFunctional(mPW91XFunctional::PW91);
|
---|
902 | funcs_[1] = new PW91CFunctional;
|
---|
903 | }
|
---|
904 | else if (!strcmp(name_,"mPWPW91")) {
|
---|
905 | init_arrays(2);
|
---|
906 | funcs_[0] = new mPW91XFunctional(mPW91XFunctional::mPW91);
|
---|
907 | funcs_[1] = new PW91CFunctional;
|
---|
908 | }
|
---|
909 | else if (!strcmp(name_,"mPW1PW91")) {
|
---|
910 | init_arrays(2);
|
---|
911 | a0_ = 0.16;
|
---|
912 | coefs_[0] = 0.84;
|
---|
913 | coefs_[1] = 1.0;
|
---|
914 | funcs_[0] = new mPW91XFunctional(mPW91XFunctional::mPW91);
|
---|
915 | funcs_[1] = new PW91CFunctional;
|
---|
916 | }
|
---|
917 | else {
|
---|
918 | ExEnv::out0() << "StdDenFunctional: bad name: " << name_ << endl;
|
---|
919 | abort();
|
---|
920 | }
|
---|
921 | }
|
---|
922 | }
|
---|
923 |
|
---|
924 | StdDenFunctional::~StdDenFunctional()
|
---|
925 | {
|
---|
926 | delete[] name_;
|
---|
927 | }
|
---|
928 |
|
---|
929 | void
|
---|
930 | StdDenFunctional::save_data_state(StateOut& s)
|
---|
931 | {
|
---|
932 | SumDenFunctional::save_data_state(s);
|
---|
933 | s.putstring(name_);
|
---|
934 | }
|
---|
935 |
|
---|
936 | void
|
---|
937 | StdDenFunctional::print(ostream& o) const
|
---|
938 | {
|
---|
939 | const char *n = name_;
|
---|
940 | if (!n) n = "Null";
|
---|
941 |
|
---|
942 | o
|
---|
943 | << indent << "Standard Density Functional: " << n << endl;
|
---|
944 | SumDenFunctional::print(o);
|
---|
945 | }
|
---|
946 |
|
---|
947 | /////////////////////////////////////////////////////////////////////////////
|
---|
948 | // LSDACFunctional: All local correlation functionals inherit from this class.
|
---|
949 | // Coded by Matt Leininger
|
---|
950 | static ClassDesc LSDACFunctional_cd(
|
---|
951 | typeid(LSDACFunctional),"LSDACFunctional",1,"public DenFunctional",
|
---|
952 | 0, 0, 0);
|
---|
953 |
|
---|
954 | LSDACFunctional::LSDACFunctional(StateIn& s):
|
---|
955 | SavableState(s),
|
---|
956 | DenFunctional(s)
|
---|
957 | {
|
---|
958 | }
|
---|
959 |
|
---|
960 | LSDACFunctional::LSDACFunctional()
|
---|
961 | {
|
---|
962 | }
|
---|
963 |
|
---|
964 | LSDACFunctional::LSDACFunctional(const Ref<KeyVal>& keyval):
|
---|
965 | DenFunctional(keyval)
|
---|
966 | {
|
---|
967 | }
|
---|
968 |
|
---|
969 | LSDACFunctional::~LSDACFunctional()
|
---|
970 | {
|
---|
971 | }
|
---|
972 |
|
---|
973 | void
|
---|
974 | LSDACFunctional::save_data_state(StateOut& s)
|
---|
975 | {
|
---|
976 | DenFunctional::save_data_state(s);
|
---|
977 | }
|
---|
978 |
|
---|
979 | void
|
---|
980 | LSDACFunctional::point(const PointInputData &id,
|
---|
981 | PointOutputData &od)
|
---|
982 | {
|
---|
983 | double junk_1, junk_2, junk_3;
|
---|
984 | point_lc(id, od, junk_1, junk_2, junk_3);
|
---|
985 | }
|
---|
986 |
|
---|
987 | /////////////////////////////////////////////////////////////////////////////
|
---|
988 | // SlaterXFunctional
|
---|
989 |
|
---|
990 | static ClassDesc SlaterXFunctional_cd(
|
---|
991 | typeid(SlaterXFunctional),"SlaterXFunctional",1,"public DenFunctional",
|
---|
992 | 0, create<SlaterXFunctional>, create<SlaterXFunctional>);
|
---|
993 |
|
---|
994 | SlaterXFunctional::SlaterXFunctional(StateIn& s):
|
---|
995 | SavableState(s),
|
---|
996 | DenFunctional(s)
|
---|
997 | {
|
---|
998 | }
|
---|
999 |
|
---|
1000 | SlaterXFunctional::SlaterXFunctional()
|
---|
1001 | {
|
---|
1002 | }
|
---|
1003 |
|
---|
1004 | SlaterXFunctional::SlaterXFunctional(const Ref<KeyVal>& keyval):
|
---|
1005 | DenFunctional(keyval)
|
---|
1006 | {
|
---|
1007 | }
|
---|
1008 |
|
---|
1009 | SlaterXFunctional::~SlaterXFunctional()
|
---|
1010 | {
|
---|
1011 | }
|
---|
1012 |
|
---|
1013 | void
|
---|
1014 | SlaterXFunctional::save_data_state(StateOut& s)
|
---|
1015 | {
|
---|
1016 | DenFunctional::save_data_state(s);
|
---|
1017 | }
|
---|
1018 |
|
---|
1019 | void
|
---|
1020 | SlaterXFunctional::point(const PointInputData &id,
|
---|
1021 | PointOutputData &od)
|
---|
1022 | {
|
---|
1023 | const double mcx2rthird = -0.9305257363491; // -1.5*(3/4pi)^1/3
|
---|
1024 | const double dmcx2rthird = -1.2407009817988; // 2*(3/4pi)^1/3
|
---|
1025 | od.zero();
|
---|
1026 |
|
---|
1027 | if (!spin_polarized_) {
|
---|
1028 | od.energy = mcx2rthird * 2.0 * id.a.rho * id.a.rho_13;
|
---|
1029 | if (compute_potential_) {
|
---|
1030 | od.df_drho_a = dmcx2rthird * id.a.rho_13;
|
---|
1031 | od.df_drho_b = od.df_drho_a;
|
---|
1032 | }
|
---|
1033 | }
|
---|
1034 | else {
|
---|
1035 | od.energy = mcx2rthird
|
---|
1036 | * (id.a.rho * id.a.rho_13 + id.b.rho * id.b.rho_13);
|
---|
1037 | if (compute_potential_) {
|
---|
1038 | od.df_drho_a = dmcx2rthird * id.a.rho_13;
|
---|
1039 | od.df_drho_b = dmcx2rthird * id.b.rho_13;
|
---|
1040 | }
|
---|
1041 | }
|
---|
1042 | }
|
---|
1043 |
|
---|
1044 | /////////////////////////////////////////////////////////////////////////////
|
---|
1045 | // PW92LCFunctional
|
---|
1046 | // Coded by Matt Leininger
|
---|
1047 | static ClassDesc PW92LCFunctional_cd(
|
---|
1048 | typeid(PW92LCFunctional),"PW92LCFunctional",1,"public LSDACFunctional",
|
---|
1049 | 0, create<PW92LCFunctional>, create<PW92LCFunctional>);
|
---|
1050 |
|
---|
1051 | PW92LCFunctional::PW92LCFunctional(StateIn& s):
|
---|
1052 | SavableState(s),
|
---|
1053 | LSDACFunctional(s)
|
---|
1054 | {
|
---|
1055 | }
|
---|
1056 |
|
---|
1057 | PW92LCFunctional::PW92LCFunctional()
|
---|
1058 | {
|
---|
1059 | }
|
---|
1060 |
|
---|
1061 | PW92LCFunctional::PW92LCFunctional(const Ref<KeyVal>& keyval):
|
---|
1062 | LSDACFunctional(keyval)
|
---|
1063 | {
|
---|
1064 | }
|
---|
1065 |
|
---|
1066 | PW92LCFunctional::~PW92LCFunctional()
|
---|
1067 | {
|
---|
1068 | }
|
---|
1069 |
|
---|
1070 | void
|
---|
1071 | PW92LCFunctional::save_data_state(StateOut& s)
|
---|
1072 | {
|
---|
1073 | LSDACFunctional::save_data_state(s);
|
---|
1074 | }
|
---|
1075 |
|
---|
1076 | double
|
---|
1077 | PW92LCFunctional::F(double x, double A, double alpha_1, double beta_1, double beta_2,
|
---|
1078 | double beta_3, double beta_4, double p)
|
---|
1079 | {
|
---|
1080 | double x2 = x*x; // r_s
|
---|
1081 | double denom = 2.*A*( beta_1 * x + beta_2 * x2 + beta_3 * x2*x + beta_4 * pow(x2,p+1.));
|
---|
1082 | double res = -2.*A*(1. + alpha_1*x2)*log(1.+ 1./denom);
|
---|
1083 |
|
---|
1084 | return res;
|
---|
1085 | }
|
---|
1086 |
|
---|
1087 | double
|
---|
1088 | PW92LCFunctional::dFdr_s(double x, double A, double alpha_1, double beta_1, double beta_2,
|
---|
1089 | double beta_3, double beta_4, double p)
|
---|
1090 | {
|
---|
1091 | double x2 = x*x; // r_s
|
---|
1092 | double Q_0 = -2.*A*(1. + alpha_1*x2);
|
---|
1093 | double Q_1 = 2.*A*(beta_1 * x + beta_2 * x2 + beta_3*x*x2 + beta_4 * pow(x2,p+1.));
|
---|
1094 | double Q_1prime = A *
|
---|
1095 | ( beta_1 * 1./x + 2.*beta_2 + 3.*beta_3*x + 2.*(p+1.)*beta_4*pow(x2,p));
|
---|
1096 | double res = -2.*A*alpha_1*log(1. + 1./Q_1) - Q_0*Q_1prime/(Q_1*Q_1 + Q_1);
|
---|
1097 |
|
---|
1098 | return res;
|
---|
1099 | }
|
---|
1100 |
|
---|
1101 | void
|
---|
1102 | PW92LCFunctional::point_lc(const PointInputData &id, PointOutputData &od,
|
---|
1103 | double &ec_local, double &decrs, double &deczeta)
|
---|
1104 | {
|
---|
1105 | od.zero();
|
---|
1106 | const double fpp0 = 4./9. * 1./(pow(2., (1./3.)) - 1.);
|
---|
1107 | const double sixth = 1./6.;
|
---|
1108 | const double four_thirds = 4./3.;
|
---|
1109 | const double one_third = 1./3.;
|
---|
1110 | const double two_thirds = 2./3.;
|
---|
1111 |
|
---|
1112 | double rho = id.a.rho + id.b.rho;
|
---|
1113 | double zeta = (id.a.rho - id.b.rho)/rho;
|
---|
1114 | double x = pow(3./(4.*M_PI*rho), sixth);
|
---|
1115 | double rs = x*x;
|
---|
1116 |
|
---|
1117 | double epc = F(x, 0.0310907, 0.21370, 7.5957, 3.5876, 1.6382, 0.49294, 1.00);
|
---|
1118 | double efc = F(x, 0.01554535, 0.20548, 14.1189, 6.1977, 3.3662, 0.62517, 1.00);
|
---|
1119 | double alphac = F(x, 0.0168869, 0.11125, 10.357, 3.6231, 0.88026, 0.49671, 1.00);
|
---|
1120 |
|
---|
1121 | double f = 9./8.*fpp0*(pow(1.+zeta, four_thirds)+pow(1.-zeta, four_thirds)-2.);
|
---|
1122 | double zeta2 = zeta*zeta;
|
---|
1123 | double zeta4 = zeta2*zeta2;
|
---|
1124 | double delta_ec = -alphac * f / fpp0 * (1. - zeta4) + (efc - epc) * f * zeta4;
|
---|
1125 | double ec = epc + delta_ec;
|
---|
1126 |
|
---|
1127 | od.energy = ec * rho;
|
---|
1128 | ec_local = ec;
|
---|
1129 |
|
---|
1130 | if (compute_potential_) {
|
---|
1131 | if (!spin_polarized_) {
|
---|
1132 | double depc_dr_s0 =
|
---|
1133 | dFdr_s(x, 0.0310907, 0.21370, 7.5957, 3.5876, 1.6382, 0.49294, 1.00);
|
---|
1134 | double dec_dr_s = depc_dr_s0;
|
---|
1135 | od.df_drho_a = od.df_drho_b = ec - (rs/3.)*dec_dr_s;
|
---|
1136 | decrs = dec_dr_s;
|
---|
1137 | deczeta = 0.;
|
---|
1138 | }
|
---|
1139 | else {
|
---|
1140 | double zeta3 = zeta2*zeta;
|
---|
1141 | double depc_dr_s0 = dFdr_s(x, 0.0310907, 0.21370, 7.5957,
|
---|
1142 | 3.5876, 1.6382, 0.49294, 1.00);
|
---|
1143 | double defc_dr_s1 = dFdr_s(x, 0.01554535, 0.20548, 14.1189,
|
---|
1144 | 6.1977, 3.3662, 0.62517, 1.00);
|
---|
1145 | double dalphac_dr_s = dFdr_s(x, 0.0168869, 0.11125, 10.357,
|
---|
1146 | 3.6231, 0.88026, 0.49671, 1.00);
|
---|
1147 | double dec_dr_s = depc_dr_s0*(1 - f*zeta4) + defc_dr_s1 * f * zeta4
|
---|
1148 | + -dalphac_dr_s * f / fpp0 * (1 - zeta4);
|
---|
1149 | double fp = two_thirds * (pow((1+zeta),one_third)
|
---|
1150 | - pow((1-zeta),one_third))/(pow(2.,one_third)-1);
|
---|
1151 | double dec_dzeta = 4.* zeta3 * f * (efc - epc - (-alphac/fpp0))
|
---|
1152 | + fp * (zeta4 * (efc - epc) + (1-zeta4)*(-alphac/fpp0));
|
---|
1153 | od.df_drho_a = ec - (rs/3.)*dec_dr_s - (zeta-1)*dec_dzeta;
|
---|
1154 | od.df_drho_b = ec - (rs/3.)*dec_dr_s - (zeta+1)*dec_dzeta;
|
---|
1155 | decrs = dec_dr_s;
|
---|
1156 | deczeta = dec_dzeta;
|
---|
1157 | }
|
---|
1158 | }
|
---|
1159 | }
|
---|
1160 |
|
---|
1161 | /////////////////////////////////////////////////////////////////////////////
|
---|
1162 | // PZ81LCFunctional
|
---|
1163 | // Coded by Matt Leininger
|
---|
1164 | // Used in P86 correlation functional
|
---|
1165 | // J. P. Perdew and A. Zunger, Phys. Rev. B, 23, 5048, 1981.
|
---|
1166 | // C. W. Murray, N. C. Handy, G. J. Laming, Mol. Phys., 78, 997, 1993.
|
---|
1167 | //
|
---|
1168 | static ClassDesc PZ81LCFunctional_cd(
|
---|
1169 | typeid(PZ81LCFunctional),"PZ81LCFunctional",1,"public LSDACFunctional",
|
---|
1170 | 0, create<PZ81LCFunctional>, create<PZ81LCFunctional>);
|
---|
1171 |
|
---|
1172 | PZ81LCFunctional::PZ81LCFunctional(StateIn& s):
|
---|
1173 | SavableState(s),
|
---|
1174 | LSDACFunctional(s)
|
---|
1175 | {
|
---|
1176 | }
|
---|
1177 |
|
---|
1178 | PZ81LCFunctional::PZ81LCFunctional()
|
---|
1179 | {
|
---|
1180 | }
|
---|
1181 |
|
---|
1182 | PZ81LCFunctional::PZ81LCFunctional(const Ref<KeyVal>& keyval):
|
---|
1183 | LSDACFunctional(keyval)
|
---|
1184 | {
|
---|
1185 | }
|
---|
1186 |
|
---|
1187 | PZ81LCFunctional::~PZ81LCFunctional()
|
---|
1188 | {
|
---|
1189 | }
|
---|
1190 |
|
---|
1191 | void
|
---|
1192 | PZ81LCFunctional::save_data_state(StateOut& s)
|
---|
1193 | {
|
---|
1194 | LSDACFunctional::save_data_state(s);
|
---|
1195 | }
|
---|
1196 |
|
---|
1197 | double
|
---|
1198 | PZ81LCFunctional::Fec_rsgt1(double rs, double beta1, double beta2, double gamma)
|
---|
1199 | {
|
---|
1200 | double sqrt_rs = sqrt(rs);
|
---|
1201 | double res = gamma / (1. + beta1*sqrt_rs + beta2*rs);
|
---|
1202 |
|
---|
1203 | return res;
|
---|
1204 | }
|
---|
1205 |
|
---|
1206 | double
|
---|
1207 | PZ81LCFunctional::dFec_rsgt1_drho(double rs, double beta1, double beta2, double gamma,
|
---|
1208 | double &dec_drs)
|
---|
1209 | {
|
---|
1210 | double ec = Fec_rsgt1(rs, beta1, beta2, gamma);
|
---|
1211 | double sqrt_rs = sqrt(rs);
|
---|
1212 | // double numer = 1.+ 7./6.*beta1*sqrt_rs + 4./3.*beta2*rs;
|
---|
1213 | double denom = 1. + beta1*sqrt_rs + beta2*rs;
|
---|
1214 | dec_drs = -ec/denom * (beta1/(2.*sqrt_rs) + beta2);
|
---|
1215 | // double res = ec * numer / denom;
|
---|
1216 | double res = (ec - rs/3.*dec_drs);
|
---|
1217 |
|
---|
1218 | return res;
|
---|
1219 | }
|
---|
1220 |
|
---|
1221 | double
|
---|
1222 | PZ81LCFunctional::Fec_rslt1(double rs, double A, double B, double C, double D)
|
---|
1223 | {
|
---|
1224 | double lnrs = log(rs);
|
---|
1225 | double res = A*lnrs + B + C*rs*lnrs + D*rs;
|
---|
1226 |
|
---|
1227 | return res;
|
---|
1228 | }
|
---|
1229 |
|
---|
1230 | double
|
---|
1231 | PZ81LCFunctional::dFec_rslt1_drho(double rs, double A, double B, double C, double D,
|
---|
1232 | double &dec_drs)
|
---|
1233 | {
|
---|
1234 | double lnrs = log(rs);
|
---|
1235 | double res = A*lnrs + B - A/3. + 2./3.*C*rs*lnrs + 1./3.*(2.*D - C)*rs;
|
---|
1236 | dec_drs = A/rs + C*lnrs + C + D;
|
---|
1237 | return res;
|
---|
1238 | }
|
---|
1239 |
|
---|
1240 |
|
---|
1241 | void
|
---|
1242 | PZ81LCFunctional::point_lc(const PointInputData &id, PointOutputData &od,
|
---|
1243 | double &ec_local, double &decrs, double &deczeta)
|
---|
1244 | {
|
---|
1245 | od.zero();
|
---|
1246 |
|
---|
1247 | const double Au = 0.0311;
|
---|
1248 | const double Ap = 0.01555;
|
---|
1249 | const double Bu = -0.048;
|
---|
1250 | const double Bp = -0.0269;
|
---|
1251 | const double Cu = 0.0020;
|
---|
1252 | const double Cp = 0.0007;
|
---|
1253 | const double Du = -0.0116;
|
---|
1254 | const double Dp = -0.0048;
|
---|
1255 | const double beta1u = 1.0529;
|
---|
1256 | const double beta1p = 1.3981;
|
---|
1257 | const double beta2u = 0.3334;
|
---|
1258 | const double beta2p = 0.2611;
|
---|
1259 | const double gammau = -0.1423;
|
---|
1260 | const double gammap = -0.0843;
|
---|
1261 | double rho = id.a.rho + id.b.rho;
|
---|
1262 | double zeta = (id.a.rho - id.b.rho)/rho;
|
---|
1263 | double rs = pow(3./(4.*M_PI*rho), (1./3.) );
|
---|
1264 | double fzeta = ( pow((1.+zeta), (4./3.)) + pow((1.-zeta), (4./3.)) - 2.)
|
---|
1265 | / ( pow(2., (4./3.)) - 2. );
|
---|
1266 |
|
---|
1267 | double euc, epc;
|
---|
1268 | if (rs >= 1.) {
|
---|
1269 | // Ceperley U
|
---|
1270 | // euc = Fec_rsgt1(rs, 1.1581, 0.3446, -0.1471);
|
---|
1271 | // Ceperley P
|
---|
1272 | // epc = Fec_rsgt1(rs, 1.2520, 0.2567, -0.0790);
|
---|
1273 | // Ceperley-Adler U
|
---|
1274 | euc = Fec_rsgt1(rs, beta1u, beta2u, gammau);
|
---|
1275 | // Ceperley-Adler P
|
---|
1276 | epc = Fec_rsgt1(rs, beta1p, beta2p, gammap);
|
---|
1277 | }
|
---|
1278 | else { // rs < 1.
|
---|
1279 | // Ceperley U with A_u and B_u
|
---|
1280 | // euc = Fec_rslt1(rs, 0.0311, -0.048, 0.0014, -0.0108);
|
---|
1281 | // Ceperley P with A_p and B_p
|
---|
1282 | // epc = Fec_rslt1(rs, 0.01555, -0.0269, 0.0001, -0.0046);
|
---|
1283 | // Ceperley-Adler U with A_u and B_u
|
---|
1284 | euc = Fec_rslt1(rs, Au, Bu, Cu, Du);
|
---|
1285 | // Ceperley-Adler P with A_p and B_p
|
---|
1286 | epc = Fec_rslt1(rs, Ap, Bp, Cp, Dp);
|
---|
1287 | }
|
---|
1288 | double ec = euc + fzeta*(epc-euc);
|
---|
1289 | ec_local = ec;
|
---|
1290 | od.energy = ec * rho;
|
---|
1291 |
|
---|
1292 | if (compute_potential_) {
|
---|
1293 | double deuc_drs = 0.;
|
---|
1294 | double depc_drs = 0.;
|
---|
1295 | double deuc_drho, depc_drho;
|
---|
1296 | if (rs > 1.) {
|
---|
1297 | // Ceperley U
|
---|
1298 | // deuc_drho = dFec_rsgt1_drho(rs, 1.1581, 0.3446, -0.1471, deuc_drs);
|
---|
1299 | // Ceperley P
|
---|
1300 | // depc_drho = dFec_rsgt1_drho(rs, 1.2520, 0.2567, -0.0790, depc_drs);
|
---|
1301 | // Ceperley-Adler U
|
---|
1302 | deuc_drho = dFec_rsgt1_drho(rs, beta1u, beta2u, gammau, deuc_drs);
|
---|
1303 | // Ceperley-Adler P
|
---|
1304 | depc_drho = dFec_rsgt1_drho(rs, beta1p, beta2p, gammap, depc_drs);
|
---|
1305 | }
|
---|
1306 | else { // rs < 1.
|
---|
1307 | // Ceperley U with A_u and B_u
|
---|
1308 | // deuc_drho = dFec_rslt1_drho(rs, 0.0311, -0.048, 0.0014, -0.0108, deuc_drs);
|
---|
1309 | // Ceperley P with A_p and B_p
|
---|
1310 | // depc_drho = dFec_rslt1_drho(rs, 0.01555, -0.0269, 0.0001, -0.0046, depc_drs);
|
---|
1311 | // Ceperley-Adler U with A_u and B_u
|
---|
1312 | deuc_drho = dFec_rslt1_drho(rs, Au, Bu, Cu, Du, deuc_drs);
|
---|
1313 | // Ceperley-Adler P with A_p and B_p
|
---|
1314 | depc_drho = dFec_rslt1_drho(rs, Ap, Bp, Cp, Dp, depc_drs);
|
---|
1315 | }
|
---|
1316 | double dfzeta_dzeta = 4./3.*( pow((1.+zeta), (1./3.)) - pow((1.-zeta), (1./3.)) )
|
---|
1317 | / (pow(2., (4./3.)) - 2.);
|
---|
1318 | decrs = deuc_drs + fzeta*(depc_drs - deuc_drs);
|
---|
1319 | deczeta = dfzeta_dzeta*(epc - euc);
|
---|
1320 | od.df_drho_a = deuc_drho + fzeta*(depc_drho - deuc_drho) +
|
---|
1321 | (epc - euc)*(1.-zeta)*dfzeta_dzeta;
|
---|
1322 | od.df_drho_b = deuc_drho + fzeta*(depc_drho - deuc_drho) +
|
---|
1323 | (epc - euc)*(-1.-zeta)*dfzeta_dzeta;
|
---|
1324 |
|
---|
1325 | }
|
---|
1326 | }
|
---|
1327 |
|
---|
1328 | /////////////////////////////////////////////////////////////////////////////
|
---|
1329 | // VWNLCFunctional
|
---|
1330 | // Coded by Matt Leininger
|
---|
1331 |
|
---|
1332 | static ClassDesc VWNLCFunctional_cd(
|
---|
1333 | typeid(VWNLCFunctional),"VWNLCFunctional",1,"public LSDACFunctional",
|
---|
1334 | 0, create<VWNLCFunctional>, create<VWNLCFunctional>);
|
---|
1335 |
|
---|
1336 | VWNLCFunctional::VWNLCFunctional(StateIn& s):
|
---|
1337 | SavableState(s), LSDACFunctional(s)
|
---|
1338 | {
|
---|
1339 | s.get(Ap_);
|
---|
1340 | s.get(Af_);
|
---|
1341 | s.get(A_alpha_);
|
---|
1342 |
|
---|
1343 | s.get(x0p_mc_);
|
---|
1344 | s.get(bp_mc_);
|
---|
1345 | s.get(cp_mc_);
|
---|
1346 | s.get(x0f_mc_);
|
---|
1347 | s.get(bf_mc_);
|
---|
1348 | s.get(cf_mc_);
|
---|
1349 | s.get(x0_alpha_mc_);
|
---|
1350 | s.get(b_alpha_mc_);
|
---|
1351 | s.get(c_alpha_mc_);
|
---|
1352 |
|
---|
1353 | s.get(x0p_rpa_);
|
---|
1354 | s.get(bp_rpa_);
|
---|
1355 | s.get(cp_rpa_);
|
---|
1356 | s.get(x0f_rpa_);
|
---|
1357 | s.get(bf_rpa_);
|
---|
1358 | s.get(cf_rpa_);
|
---|
1359 | s.get(x0_alpha_rpa_);
|
---|
1360 | s.get(b_alpha_rpa_);
|
---|
1361 | s.get(c_alpha_rpa_);
|
---|
1362 | }
|
---|
1363 |
|
---|
1364 | VWNLCFunctional::VWNLCFunctional()
|
---|
1365 | {
|
---|
1366 | init_constants();
|
---|
1367 | }
|
---|
1368 |
|
---|
1369 | VWNLCFunctional::~VWNLCFunctional()
|
---|
1370 | {
|
---|
1371 | }
|
---|
1372 |
|
---|
1373 | VWNLCFunctional::VWNLCFunctional(const Ref<KeyVal>& keyval):
|
---|
1374 | LSDACFunctional(keyval)
|
---|
1375 | {
|
---|
1376 | init_constants();
|
---|
1377 | Ap_ = keyval->doublevalue("Ap", KeyValValuedouble(Ap_));
|
---|
1378 | Af_ = keyval->doublevalue("Af", KeyValValuedouble(Af_));
|
---|
1379 | A_alpha_ = keyval->doublevalue("A_alpha", KeyValValuedouble(A_alpha_));
|
---|
1380 |
|
---|
1381 | x0p_mc_ = keyval->doublevalue("x0p_mc", KeyValValuedouble(x0p_mc_));
|
---|
1382 | bp_mc_ = keyval->doublevalue("bp_mc", KeyValValuedouble(bp_mc_));
|
---|
1383 | cp_mc_ = keyval->doublevalue("cp_mc", KeyValValuedouble(cp_mc_));
|
---|
1384 | x0f_mc_ = keyval->doublevalue("x0f_mc", KeyValValuedouble(x0f_mc_));
|
---|
1385 | bf_mc_ = keyval->doublevalue("bf_mc", KeyValValuedouble(bf_mc_));
|
---|
1386 | cf_mc_ = keyval->doublevalue("cf_mc", KeyValValuedouble(cf_mc_));
|
---|
1387 | x0_alpha_mc_ = keyval->doublevalue("x0_alpha_mc",
|
---|
1388 | KeyValValuedouble(x0_alpha_mc_));
|
---|
1389 | b_alpha_mc_ = keyval->doublevalue("b_alpha_mc",
|
---|
1390 | KeyValValuedouble(b_alpha_mc_));
|
---|
1391 | c_alpha_mc_ = keyval->doublevalue("c_alpha_mc",
|
---|
1392 | KeyValValuedouble(c_alpha_mc_));
|
---|
1393 |
|
---|
1394 | x0p_rpa_ = keyval->doublevalue("x0p_rpa", KeyValValuedouble(x0p_rpa_));
|
---|
1395 | bp_rpa_ = keyval->doublevalue("bp_rpa", KeyValValuedouble(bp_rpa_));
|
---|
1396 | cp_rpa_ = keyval->doublevalue("cp_rpa", KeyValValuedouble(cp_rpa_));
|
---|
1397 | x0f_rpa_ = keyval->doublevalue("x0f_rpa", KeyValValuedouble(x0f_rpa_));
|
---|
1398 | bf_rpa_ = keyval->doublevalue("bf_rpa", KeyValValuedouble(bf_rpa_));
|
---|
1399 | cf_rpa_ = keyval->doublevalue("cf_rpa", KeyValValuedouble(cf_rpa_));
|
---|
1400 | x0_alpha_rpa_ = keyval->doublevalue("x0_alpha_rpa",
|
---|
1401 | KeyValValuedouble(x0_alpha_rpa_));
|
---|
1402 | b_alpha_rpa_ = keyval->doublevalue("b_alpha_rpa",
|
---|
1403 | KeyValValuedouble(b_alpha_rpa_));
|
---|
1404 | c_alpha_rpa_ = keyval->doublevalue("c_alpha_rpa",
|
---|
1405 | KeyValValuedouble(c_alpha_rpa_));
|
---|
1406 | }
|
---|
1407 |
|
---|
1408 | void
|
---|
1409 | VWNLCFunctional::save_data_state(StateOut& s)
|
---|
1410 | {
|
---|
1411 | LSDACFunctional::save_data_state(s);
|
---|
1412 |
|
---|
1413 | s.put(Ap_);
|
---|
1414 | s.put(Af_);
|
---|
1415 | s.put(A_alpha_);
|
---|
1416 |
|
---|
1417 | s.put(x0p_mc_);
|
---|
1418 | s.put(bp_mc_);
|
---|
1419 | s.put(cp_mc_);
|
---|
1420 | s.put(x0f_mc_);
|
---|
1421 | s.put(bf_mc_);
|
---|
1422 | s.put(cf_mc_);
|
---|
1423 | s.put(x0_alpha_mc_);
|
---|
1424 | s.put(b_alpha_mc_);
|
---|
1425 | s.put(c_alpha_mc_);
|
---|
1426 |
|
---|
1427 | s.put(x0p_rpa_);
|
---|
1428 | s.put(bp_rpa_);
|
---|
1429 | s.put(cp_rpa_);
|
---|
1430 | s.put(x0f_rpa_);
|
---|
1431 | s.put(bf_rpa_);
|
---|
1432 | s.put(cf_rpa_);
|
---|
1433 | s.put(x0_alpha_rpa_);
|
---|
1434 | s.put(b_alpha_rpa_);
|
---|
1435 | s.put(c_alpha_rpa_);
|
---|
1436 | }
|
---|
1437 |
|
---|
1438 | void
|
---|
1439 | VWNLCFunctional::init_constants()
|
---|
1440 | {
|
---|
1441 | Ap_ = 0.0310907;
|
---|
1442 | Af_ = 0.01554535;
|
---|
1443 | A_alpha_ = -1./(6.*M_PI*M_PI);
|
---|
1444 |
|
---|
1445 | x0p_mc_ = -0.10498;
|
---|
1446 | bp_mc_ = 3.72744;
|
---|
1447 | cp_mc_ = 12.9352;
|
---|
1448 | x0f_mc_ = -0.32500;
|
---|
1449 | bf_mc_ = 7.06042;
|
---|
1450 | cf_mc_ = 18.0578;
|
---|
1451 | x0_alpha_mc_ = -0.00475840;
|
---|
1452 | b_alpha_mc_ = 1.13107;
|
---|
1453 | c_alpha_mc_ = 13.0045;
|
---|
1454 |
|
---|
1455 | x0p_rpa_ = -0.409286;
|
---|
1456 | bp_rpa_ = 13.0720;
|
---|
1457 | cp_rpa_ = 42.7198;
|
---|
1458 | x0f_rpa_ = -0.743294;
|
---|
1459 | bf_rpa_ = 20.1231;
|
---|
1460 | cf_rpa_ = 101.578;
|
---|
1461 | x0_alpha_rpa_ = -0.228344;
|
---|
1462 | b_alpha_rpa_ = 1.06835;
|
---|
1463 | c_alpha_rpa_ = 11.4813;
|
---|
1464 | }
|
---|
1465 |
|
---|
1466 | double
|
---|
1467 | VWNLCFunctional::F(double x, double A, double x0, double b, double c)
|
---|
1468 | {
|
---|
1469 | double x2 = x*x;
|
---|
1470 | double x02 = x0*x0;
|
---|
1471 | double Xx = x2 + b*x + c;
|
---|
1472 | double Xx0 = x02 + b*x0 + c;
|
---|
1473 | double Q = sqrt(4.*c-b*b);
|
---|
1474 | double res
|
---|
1475 | = A * ( log(x2/Xx)
|
---|
1476 | + 2.*b/Q * atan(Q/(2.*x+b))
|
---|
1477 | - b*x0/Xx0 * ( log((x-x0)*(x-x0)/Xx)
|
---|
1478 | + 2.*(b+2.*x0)/Q * atan(Q/(2.*x+b)) ) );
|
---|
1479 | return res;
|
---|
1480 | }
|
---|
1481 |
|
---|
1482 | double
|
---|
1483 | VWNLCFunctional::dFdr_s(double x, double A, double x0, double b, double c)
|
---|
1484 | {
|
---|
1485 | double x2 = x*x;
|
---|
1486 | double x02 = x0*x0;
|
---|
1487 | double Xx = x2 + b*x +c;
|
---|
1488 | double Xx0 = x02 + b*x0 + c;
|
---|
1489 | double Q = sqrt(4.*c-b*b);
|
---|
1490 | double res
|
---|
1491 | = A * ( 1./x2 - 1./Xx - b/(2.*Xx*x)
|
---|
1492 | + ((x0*(2.*x0+b))/Xx0 - 1) * (2.*b)/(x*(Q*Q+(2.*x+b)*(2.*x+b)))
|
---|
1493 | - (b*x0)/(x*(x-x0)*Xx0) + (b*x0*(1+(b/(2.*x))))/(Xx0*Xx) );
|
---|
1494 | return res;
|
---|
1495 | }
|
---|
1496 |
|
---|
1497 | void
|
---|
1498 | VWNLCFunctional::point_lc(const PointInputData &id, PointOutputData &od,
|
---|
1499 | double &ec_local, double &decrs, double &deczeta)
|
---|
1500 | {
|
---|
1501 | }
|
---|
1502 |
|
---|
1503 | /////////////////////////////////////////////////////////////////////////////
|
---|
1504 | // VWN1LCFunctional
|
---|
1505 | // Coded by Matt Leininger
|
---|
1506 |
|
---|
1507 | static ClassDesc VWN1LCFunctional_cd(
|
---|
1508 | typeid(VWN1LCFunctional),"VWN1LCFunctional",1,"public LSDACFunctional",
|
---|
1509 | 0, create<VWN1LCFunctional>, create<VWN1LCFunctional>);
|
---|
1510 |
|
---|
1511 | VWN1LCFunctional::VWN1LCFunctional(StateIn& s):
|
---|
1512 | SavableState(s),
|
---|
1513 | VWNLCFunctional(s)
|
---|
1514 | {
|
---|
1515 | s.get(x0p_);
|
---|
1516 | s.get(bp_);
|
---|
1517 | s.get(cp_);
|
---|
1518 | s.get(x0f_);
|
---|
1519 | s.get(bf_);
|
---|
1520 | s.get(cf_);
|
---|
1521 | }
|
---|
1522 |
|
---|
1523 | VWN1LCFunctional::VWN1LCFunctional()
|
---|
1524 | {
|
---|
1525 | x0p_ = x0p_mc_;
|
---|
1526 | bp_ = bp_mc_;
|
---|
1527 | cp_ = cp_mc_;
|
---|
1528 | x0f_ = x0f_mc_;
|
---|
1529 | bf_ = bf_mc_;
|
---|
1530 | cf_ = cf_mc_;
|
---|
1531 | }
|
---|
1532 |
|
---|
1533 | VWN1LCFunctional::VWN1LCFunctional(int use_rpa)
|
---|
1534 | {
|
---|
1535 | if (use_rpa) {
|
---|
1536 | x0p_ = x0p_rpa_;
|
---|
1537 | bp_ = bp_rpa_;
|
---|
1538 | cp_ = cp_rpa_;
|
---|
1539 | x0f_ = x0f_rpa_;
|
---|
1540 | bf_ = bf_rpa_;
|
---|
1541 | cf_ = cf_rpa_;
|
---|
1542 | }
|
---|
1543 | else {
|
---|
1544 | x0p_ = x0p_mc_;
|
---|
1545 | bp_ = bp_mc_;
|
---|
1546 | cp_ = cp_mc_;
|
---|
1547 | x0f_ = x0f_mc_;
|
---|
1548 | bf_ = bf_mc_;
|
---|
1549 | cf_ = cf_mc_;
|
---|
1550 | }
|
---|
1551 |
|
---|
1552 | }
|
---|
1553 |
|
---|
1554 | VWN1LCFunctional::VWN1LCFunctional(const Ref<KeyVal>& keyval):
|
---|
1555 | VWNLCFunctional(keyval)
|
---|
1556 | {
|
---|
1557 | int vwn1rpa = keyval->booleanvalue("rpa", KeyValValueboolean(0));
|
---|
1558 | if (vwn1rpa) {
|
---|
1559 | x0p_ = x0p_rpa_;
|
---|
1560 | bp_ = bp_rpa_;
|
---|
1561 | cp_ = cp_rpa_;
|
---|
1562 | x0f_ = x0f_rpa_;
|
---|
1563 | bf_ = bf_rpa_;
|
---|
1564 | cf_ = cf_rpa_;
|
---|
1565 | }
|
---|
1566 | else {
|
---|
1567 | x0p_ = x0p_mc_;
|
---|
1568 | bp_ = bp_mc_;
|
---|
1569 | cp_ = cp_mc_;
|
---|
1570 | x0f_ = x0f_mc_;
|
---|
1571 | bf_ = bf_mc_;
|
---|
1572 | cf_ = cf_mc_;
|
---|
1573 | }
|
---|
1574 |
|
---|
1575 | x0p_ = keyval->doublevalue("x0p", KeyValValuedouble(x0p_));
|
---|
1576 | bp_ = keyval->doublevalue("bp", KeyValValuedouble(bp_));
|
---|
1577 | cp_ = keyval->doublevalue("cp", KeyValValuedouble(cp_));
|
---|
1578 | x0f_ = keyval->doublevalue("x0f", KeyValValuedouble(x0f_));
|
---|
1579 | bf_ = keyval->doublevalue("bf", KeyValValuedouble(bf_));
|
---|
1580 | cf_ = keyval->doublevalue("cf", KeyValValuedouble(cf_));
|
---|
1581 | }
|
---|
1582 |
|
---|
1583 | VWN1LCFunctional::~VWN1LCFunctional()
|
---|
1584 | {
|
---|
1585 | }
|
---|
1586 |
|
---|
1587 | void
|
---|
1588 | VWN1LCFunctional::save_data_state(StateOut& s)
|
---|
1589 | {
|
---|
1590 | VWNLCFunctional::save_data_state(s);
|
---|
1591 | s.put(x0p_);
|
---|
1592 | s.put(bp_);
|
---|
1593 | s.put(cp_);
|
---|
1594 | s.put(x0f_);
|
---|
1595 | s.put(bf_);
|
---|
1596 | s.put(cf_);
|
---|
1597 | }
|
---|
1598 |
|
---|
1599 | // Based on the VWN1 functional in Vosko, Wilk, and Nusair, Can. J. Phys.
|
---|
1600 | // 58, 1200, (1980).
|
---|
1601 | void
|
---|
1602 | VWN1LCFunctional::point_lc(const PointInputData &id, PointOutputData &od,
|
---|
1603 | double &ec_local, double &decrs, double &deczeta)
|
---|
1604 | {
|
---|
1605 | od.zero();
|
---|
1606 |
|
---|
1607 | double four_thirds = 4./3.;
|
---|
1608 | double one_third = 1./3.;
|
---|
1609 | double two_thirds = 2./3.;
|
---|
1610 | double sixth = 1./6.;
|
---|
1611 | double rho = id.a.rho + id.b.rho;
|
---|
1612 | double zeta = (id.a.rho - id.b.rho)/rho;
|
---|
1613 | double x = pow(3./(4.*M_PI*rho), sixth);
|
---|
1614 | double rs = x*x;
|
---|
1615 |
|
---|
1616 | double epc = F(x, Ap_, x0p_, bp_, cp_);
|
---|
1617 | double efc = F(x, Af_, x0f_, bf_, cf_);
|
---|
1618 |
|
---|
1619 | const double fpp0 = 4./9. * 1./(pow(2., (1./3.)) - 1.);
|
---|
1620 | double f = 9./8.*fpp0*(pow(1.+zeta, four_thirds)+pow(1.-zeta, four_thirds)-2.);
|
---|
1621 | double delta_ec = f * (efc - epc);
|
---|
1622 | double ec = epc + delta_ec;
|
---|
1623 |
|
---|
1624 | ec_local = ec;
|
---|
1625 | od.energy = ec * rho;
|
---|
1626 |
|
---|
1627 | if (compute_potential_) {
|
---|
1628 | if (!spin_polarized_) {
|
---|
1629 | double depc_dr_s0 = dFdr_s(x, Ap_, x0p_, bp_, cp_);
|
---|
1630 | double dec_dr_s = depc_dr_s0;
|
---|
1631 | od.df_drho_a = od.df_drho_b = ec - (rs/3.)*dec_dr_s;
|
---|
1632 | decrs = dec_dr_s;
|
---|
1633 | deczeta = 0.;
|
---|
1634 | }
|
---|
1635 | else {
|
---|
1636 | double depc_dr_s0 = dFdr_s(x, Ap_, x0p_, bp_, cp_);
|
---|
1637 | double defc_dr_s1 = dFdr_s(x, Af_, x0f_, bf_, cf_);
|
---|
1638 | double dec_dr_s = depc_dr_s0 + f * (defc_dr_s1 - depc_dr_s0);
|
---|
1639 | double fp = two_thirds * (pow((1+zeta),one_third)
|
---|
1640 | - pow((1-zeta),one_third))/(pow(2.,one_third)-1);
|
---|
1641 | // double dec_dzeta = depc_dr_s0 + fp * (efc - epc) + f * (defc_dr_s1 - depc_dr_s0);
|
---|
1642 | double dec_dzeta = fp * (efc - epc);
|
---|
1643 | od.df_drho_a = ec - (rs/3.)*dec_dr_s - (zeta-1)*dec_dzeta;
|
---|
1644 | od.df_drho_b = ec - (rs/3.)*dec_dr_s - (zeta+1)*dec_dzeta;
|
---|
1645 | decrs = dec_dr_s;
|
---|
1646 | deczeta = dec_dzeta;
|
---|
1647 | }
|
---|
1648 | }
|
---|
1649 | }
|
---|
1650 |
|
---|
1651 | /////////////////////////////////////////////////////////////////////////////
|
---|
1652 | // VWN2LCFunctional
|
---|
1653 | // Coded by Matt Leininger
|
---|
1654 | static ClassDesc VWN2LCFunctional_cd(
|
---|
1655 | typeid(VWN2LCFunctional),"VWN2LCFunctional",1,"public VWNLCFunctional",
|
---|
1656 | 0, create<VWN2LCFunctional>, create<VWN2LCFunctional>);
|
---|
1657 |
|
---|
1658 | VWN2LCFunctional::VWN2LCFunctional(StateIn& s):
|
---|
1659 | SavableState(s),
|
---|
1660 | VWNLCFunctional(s)
|
---|
1661 | {
|
---|
1662 | }
|
---|
1663 |
|
---|
1664 | VWN2LCFunctional::VWN2LCFunctional()
|
---|
1665 | {
|
---|
1666 | }
|
---|
1667 |
|
---|
1668 | VWN2LCFunctional::VWN2LCFunctional(const Ref<KeyVal>& keyval):
|
---|
1669 | VWNLCFunctional(keyval)
|
---|
1670 | {
|
---|
1671 | }
|
---|
1672 |
|
---|
1673 | VWN2LCFunctional::~VWN2LCFunctional()
|
---|
1674 | {
|
---|
1675 | }
|
---|
1676 |
|
---|
1677 | void
|
---|
1678 | VWN2LCFunctional::save_data_state(StateOut& s)
|
---|
1679 | {
|
---|
1680 | VWNLCFunctional::save_data_state(s);
|
---|
1681 | }
|
---|
1682 |
|
---|
1683 | void
|
---|
1684 | VWN2LCFunctional::point_lc(const PointInputData &id, PointOutputData &od,
|
---|
1685 | double &ec_local, double &decrs, double &deczeta)
|
---|
1686 | {
|
---|
1687 | od.zero();
|
---|
1688 | const double fpp0 = 4./9. * 1./(pow(2., (1./3.)) - 1.);
|
---|
1689 | const double sixth = 1./6.;
|
---|
1690 | const double four_thirds = 4./3.;
|
---|
1691 | const double one_third = 1./3.;
|
---|
1692 | const double two_thirds = 2./3.;
|
---|
1693 |
|
---|
1694 | double rho = id.a.rho + id.b.rho;
|
---|
1695 | double zeta = (id.a.rho - id.b.rho)/rho;
|
---|
1696 | double x = pow(3./(4.*M_PI*rho), sixth);
|
---|
1697 | double rs = x*x;
|
---|
1698 |
|
---|
1699 | // Monte Carlo fitting parameters
|
---|
1700 | double epc_mc = F(x, Ap_, x0p_mc_, bp_mc_, cp_mc_);
|
---|
1701 | double efc_mc = F(x, Af_, x0f_mc_, bf_mc_, cf_mc_);
|
---|
1702 | //double alphac_mc = F(x, A_alpha_, x0_alpha_mc_, b_alpha_mc_, c_alpha_mc_);
|
---|
1703 | // RPA fitting parameters
|
---|
1704 | double epc_rpa = F(x, Ap_, x0p_rpa_, bp_rpa_, cp_rpa_);
|
---|
1705 | double efc_rpa = F(x, Af_, x0f_rpa_, bf_rpa_, cf_rpa_);
|
---|
1706 | double alphac_rpa = F(x, A_alpha_, x0_alpha_rpa_, b_alpha_rpa_, c_alpha_rpa_);
|
---|
1707 |
|
---|
1708 | double f = 9./8.*fpp0*(pow(1.+zeta, four_thirds)+pow(1.-zeta, four_thirds)-2.);
|
---|
1709 | double zeta2 = zeta*zeta;
|
---|
1710 | double zeta4 = zeta2*zeta2;
|
---|
1711 | double delta_e_rpa = efc_rpa - epc_rpa;
|
---|
1712 | double delta_e_mc = efc_mc - epc_mc;
|
---|
1713 | double beta = (fpp0*delta_e_rpa / alphac_rpa) - 1.;
|
---|
1714 | double delta_erpa_rszeta = alphac_rpa * f / fpp0 * (1. + beta * zeta4);
|
---|
1715 | double delta_ec = delta_erpa_rszeta + f*(delta_e_mc - delta_e_rpa);
|
---|
1716 | double ec = epc_mc + delta_ec;
|
---|
1717 |
|
---|
1718 | od.energy = ec * rho;
|
---|
1719 | ec_local = ec;
|
---|
1720 |
|
---|
1721 | if (compute_potential_) {
|
---|
1722 | //double zeta3 = zeta2*zeta;
|
---|
1723 | // Monte Carlo fitting parameters
|
---|
1724 | double depc_dr_s0_mc = dFdr_s(x, Ap_,x0p_mc_, bp_mc_, cp_mc_);
|
---|
1725 | double defc_dr_s1_mc = dFdr_s(x, Af_, x0f_mc_, bf_mc_, cf_mc_);
|
---|
1726 | double dalphac_dr_s_mc = dFdr_s(x, A_alpha_, x0_alpha_mc_, b_alpha_mc_, c_alpha_mc_);
|
---|
1727 | // RPA fitting parameters
|
---|
1728 | double depc_dr_s0_rpa = dFdr_s(x, Ap_, x0p_rpa_, bp_rpa_, cp_rpa_);
|
---|
1729 | double defc_dr_s1_rpa = dFdr_s(x, Af_, x0f_rpa_, bf_rpa_, cf_rpa_);
|
---|
1730 | double dalphac_dr_s_rpa = dFdr_s(x, A_alpha_, x0_alpha_rpa_, b_alpha_rpa_, c_alpha_rpa_);
|
---|
1731 |
|
---|
1732 | double fp = two_thirds * (pow((1+zeta),one_third)
|
---|
1733 | - pow((1-zeta),one_third))/(pow(2.,one_third)-1);
|
---|
1734 |
|
---|
1735 | // RPA fitting parameters
|
---|
1736 | double ddelta_e_rpa = defc_dr_s1_rpa - depc_dr_s0_rpa;
|
---|
1737 | double ddeltae_rpa_dr_s = f / fpp0 * (1 - zeta4)* dalphac_dr_s_rpa
|
---|
1738 | + f * zeta4 * ddelta_e_rpa;
|
---|
1739 | double ddeltae_rpa_dzeta = alphac_rpa / fpp0 *
|
---|
1740 | ( fp * (1.-zeta4) - 4.* f * zeta*zeta2)
|
---|
1741 | + delta_e_rpa * ( fp*zeta4 + 4.*f*zeta*zeta2);
|
---|
1742 |
|
---|
1743 | // Monte Carlo fitting parameters
|
---|
1744 | double ddelta_e_mc = defc_dr_s1_rpa - depc_dr_s0_mc;
|
---|
1745 | // double dec_dr_s_mc = depc_dr_s0*(1 - f*zeta4) + defc_dr_s1 * f * zeta4
|
---|
1746 | // + dalphac_dr_s * f / fpp0 * (1 - zeta4);
|
---|
1747 | // double dec_dzeta_mc = 4.* zeta3 * f * (efc_mc - epc_mc - (alphac_mc/fpp0))
|
---|
1748 | // + fp * (zeta4 * (efc_mc - epc_mc) + (1-zeta4)*(alphac_mc/fpp0));
|
---|
1749 |
|
---|
1750 | double dec_dr_s = depc_dr_s0_mc + ddeltae_rpa_dr_s + f * (ddelta_e_mc - ddelta_e_rpa);
|
---|
1751 |
|
---|
1752 | double dec_dzeta = ddeltae_rpa_dzeta + fp * (delta_e_mc - delta_e_rpa);
|
---|
1753 |
|
---|
1754 | od.df_drho_a = ec - (rs/3.)*dec_dr_s - (zeta-1)*dec_dzeta;
|
---|
1755 | od.df_drho_b = ec - (rs/3.)*dec_dr_s - (zeta+1)*dec_dzeta;
|
---|
1756 | decrs = dec_dr_s;
|
---|
1757 | deczeta = dec_dzeta;
|
---|
1758 | }
|
---|
1759 | }
|
---|
1760 |
|
---|
1761 | /////////////////////////////////////////////////////////////////////////////
|
---|
1762 | // VWN3LCFunctional
|
---|
1763 | // Coded by Matt Leininger
|
---|
1764 | static ClassDesc VWN3LCFunctional_cd(
|
---|
1765 | typeid(VWN3LCFunctional),"VWN3LCFunctional",1,"public VWNLCFunctional",
|
---|
1766 | 0, create<VWN3LCFunctional>, create<VWN3LCFunctional>);
|
---|
1767 |
|
---|
1768 | VWN3LCFunctional::VWN3LCFunctional(StateIn& s):
|
---|
1769 | SavableState(s),
|
---|
1770 | VWNLCFunctional(s)
|
---|
1771 | {
|
---|
1772 | s.get(monte_carlo_prefactor_);
|
---|
1773 | s.get(monte_carlo_e0_);
|
---|
1774 | }
|
---|
1775 |
|
---|
1776 | VWN3LCFunctional::VWN3LCFunctional(int mcp, int mce0)
|
---|
1777 | {
|
---|
1778 | monte_carlo_prefactor_ = mcp;
|
---|
1779 | monte_carlo_e0_ = mce0;
|
---|
1780 | }
|
---|
1781 |
|
---|
1782 | VWN3LCFunctional::VWN3LCFunctional(const Ref<KeyVal>& keyval):
|
---|
1783 | VWNLCFunctional(keyval)
|
---|
1784 | {
|
---|
1785 | monte_carlo_prefactor_ = keyval->booleanvalue("monte_carlo_prefactor",
|
---|
1786 | KeyValValueboolean(1));
|
---|
1787 | monte_carlo_e0_ = keyval->booleanvalue("monte_carlo_e0",
|
---|
1788 | KeyValValueboolean(1));
|
---|
1789 | }
|
---|
1790 |
|
---|
1791 | VWN3LCFunctional::~VWN3LCFunctional()
|
---|
1792 | {
|
---|
1793 | }
|
---|
1794 |
|
---|
1795 | void
|
---|
1796 | VWN3LCFunctional::save_data_state(StateOut& s)
|
---|
1797 | {
|
---|
1798 | VWNLCFunctional::save_data_state(s);
|
---|
1799 | s.put(monte_carlo_prefactor_);
|
---|
1800 | s.put(monte_carlo_e0_);
|
---|
1801 | }
|
---|
1802 |
|
---|
1803 | // based on the equations given on a NIST WWW site
|
---|
1804 | void
|
---|
1805 | VWN3LCFunctional::point_lc(const PointInputData &id, PointOutputData &od,
|
---|
1806 | double &ec_local, double &decrs, double &deczeta)
|
---|
1807 | {
|
---|
1808 | od.zero();
|
---|
1809 | const double fpp0 = 4./9. * 1./(pow(2., (1./3.)) - 1.);
|
---|
1810 | const double sixth = 1./6.;
|
---|
1811 | const double four_thirds = 4./3.;
|
---|
1812 | const double one_third = 1./3.;
|
---|
1813 | const double two_thirds = 2./3.;
|
---|
1814 |
|
---|
1815 | double rho = id.a.rho + id.b.rho;
|
---|
1816 | double zeta = (id.a.rho - id.b.rho)/rho;
|
---|
1817 | double x = pow(3./(4.*M_PI*rho), sixth);
|
---|
1818 | double rs = x*x;
|
---|
1819 |
|
---|
1820 | // Monte Carlo fitting parameters
|
---|
1821 | double epc0_mc = F(x, Ap_, x0p_mc_, bp_mc_, cp_mc_);
|
---|
1822 | double efc1_mc = F(x, Af_, x0f_mc_, bf_mc_, cf_mc_);
|
---|
1823 | // double alphac_mc = F(x, A_alpha_, x0_alpha_mc_, b_alpha_mc_, c_alpha_mc_);
|
---|
1824 | // RPA fitting parameters
|
---|
1825 | double epc0_rpa = F(x, Ap_, x0p_rpa_, bp_rpa_, cp_rpa_);
|
---|
1826 | double efc1_rpa = F(x, Af_, x0f_rpa_, bf_rpa_, cf_rpa_);
|
---|
1827 | double alphac_rpa = F(x, A_alpha_, x0_alpha_rpa_, b_alpha_rpa_, c_alpha_rpa_);
|
---|
1828 |
|
---|
1829 | double f = 9./8.*fpp0*( pow(1.+zeta, four_thirds) + pow(1.-zeta, four_thirds) - 2.);
|
---|
1830 | double zeta2 = zeta*zeta;
|
---|
1831 | double zeta4 = zeta2*zeta2;
|
---|
1832 | double delta_e_rpa = efc1_rpa - epc0_rpa;
|
---|
1833 | double delta_e_mc = efc1_mc - epc0_mc;
|
---|
1834 | double delta_erpa_rszeta = alphac_rpa * f / fpp0 * (1.-zeta4) + f * zeta4 * delta_e_rpa;
|
---|
1835 | double delta_ec;
|
---|
1836 | if (!monte_carlo_prefactor_) delta_ec = delta_erpa_rszeta;
|
---|
1837 | else delta_ec = delta_e_mc/delta_e_rpa * delta_erpa_rszeta;
|
---|
1838 |
|
---|
1839 | double ec;
|
---|
1840 | if (monte_carlo_e0_) ec = epc0_mc;
|
---|
1841 | else ec = epc0_rpa;
|
---|
1842 | ec += delta_ec;
|
---|
1843 |
|
---|
1844 | od.energy = ec * rho;
|
---|
1845 | ec_local = ec;
|
---|
1846 |
|
---|
1847 | if (compute_potential_) {
|
---|
1848 | // Monte Carlo fitting parameters
|
---|
1849 | double depc_dr_s0_mc = dFdr_s(x, Ap_, x0p_mc_, bp_mc_, cp_mc_);
|
---|
1850 | double defc_dr_s1_mc = dFdr_s(x, Af_, x0f_mc_, bf_mc_, cf_mc_);
|
---|
1851 | // double dalphac_dr_s_mc = dFdr_s(x, A_alpha_, x0_alpha_mc_, b_alpha_mc_, c_alpha_mc_);
|
---|
1852 | // RPA fitting parameters
|
---|
1853 | double depc_dr_s0_rpa = dFdr_s(x, Ap_, x0p_rpa_, bp_rpa_, cp_rpa_);
|
---|
1854 | double defc_dr_s1_rpa = dFdr_s(x, Af_, x0f_rpa_, bf_rpa_, cf_rpa_);
|
---|
1855 | double dalphac_dr_s_rpa = dFdr_s(x, A_alpha_, x0_alpha_rpa_, b_alpha_rpa_, c_alpha_rpa_);
|
---|
1856 |
|
---|
1857 | double fp = two_thirds * (pow((1+zeta),one_third)
|
---|
1858 | - pow((1-zeta),one_third))/(pow(2.,one_third)-1);
|
---|
1859 |
|
---|
1860 | // RPA fitting parameters
|
---|
1861 | double ddelta_e_rpa = defc_dr_s1_rpa - depc_dr_s0_rpa;
|
---|
1862 | double ddeltae_rpa_dr_s = f / fpp0 * (1 - zeta4)* dalphac_dr_s_rpa
|
---|
1863 | + f * zeta4 * ddelta_e_rpa;
|
---|
1864 | double ddeltae_rpa_dzeta = alphac_rpa / fpp0 *
|
---|
1865 | ( fp * (1.-zeta4) - 4.* f * zeta*zeta2)
|
---|
1866 | + delta_e_rpa * ( fp*zeta4 + 4.*f*zeta*zeta2 );
|
---|
1867 |
|
---|
1868 | // Monte Carlo fitting parameters
|
---|
1869 | double ddelta_e_mc = defc_dr_s1_mc - depc_dr_s0_mc;
|
---|
1870 |
|
---|
1871 | double dec_dzeta, dec_dr_s;
|
---|
1872 | if (!monte_carlo_prefactor_) {
|
---|
1873 | dec_dzeta = ddeltae_rpa_dzeta;
|
---|
1874 | if (monte_carlo_e0_) dec_dr_s = depc_dr_s0_mc;
|
---|
1875 | else dec_dr_s = depc_dr_s0_rpa;
|
---|
1876 | dec_dr_s += ddeltae_rpa_dr_s;
|
---|
1877 | }
|
---|
1878 | else {
|
---|
1879 | dec_dzeta = delta_e_mc / delta_e_rpa * ddeltae_rpa_dzeta;
|
---|
1880 | if (monte_carlo_e0_) dec_dr_s = depc_dr_s0_mc;
|
---|
1881 | else dec_dr_s = depc_dr_s0_rpa;
|
---|
1882 | dec_dr_s += delta_erpa_rszeta/delta_e_rpa * ddelta_e_mc
|
---|
1883 | + delta_e_mc/delta_e_rpa * ddeltae_rpa_dr_s
|
---|
1884 | - delta_erpa_rszeta*delta_e_mc/(delta_e_rpa*delta_e_rpa)
|
---|
1885 | * ddelta_e_rpa;
|
---|
1886 | }
|
---|
1887 |
|
---|
1888 | od.df_drho_a = ec - (rs/3.)*dec_dr_s - (zeta-1)*dec_dzeta;
|
---|
1889 | od.df_drho_b = ec - (rs/3.)*dec_dr_s - (zeta+1)*dec_dzeta;
|
---|
1890 | decrs = dec_dr_s;
|
---|
1891 | deczeta = dec_dzeta;
|
---|
1892 | }
|
---|
1893 | }
|
---|
1894 |
|
---|
1895 | /////////////////////////////////////////////////////////////////////////////
|
---|
1896 | // VWN4LCFunctional
|
---|
1897 | // Coded by Matt Leininger
|
---|
1898 | static ClassDesc VWN4LCFunctional_cd(
|
---|
1899 | typeid(VWN4LCFunctional),"VWN4LCFunctional",1,"public VWNLCFunctional",
|
---|
1900 | 0, create<VWN4LCFunctional>, create<VWN4LCFunctional>);
|
---|
1901 |
|
---|
1902 | VWN4LCFunctional::VWN4LCFunctional(StateIn& s):
|
---|
1903 | SavableState(s),
|
---|
1904 | VWNLCFunctional(s)
|
---|
1905 | {
|
---|
1906 | s.get(monte_carlo_prefactor_);
|
---|
1907 | }
|
---|
1908 |
|
---|
1909 | VWN4LCFunctional::VWN4LCFunctional()
|
---|
1910 | {
|
---|
1911 | monte_carlo_prefactor_ = 0;
|
---|
1912 | }
|
---|
1913 |
|
---|
1914 | VWN4LCFunctional::VWN4LCFunctional(const Ref<KeyVal>& keyval):
|
---|
1915 | VWNLCFunctional(keyval)
|
---|
1916 | {
|
---|
1917 | monte_carlo_prefactor_ = keyval->booleanvalue("monte_carlo_prefactor",
|
---|
1918 | KeyValValueboolean(0));
|
---|
1919 | }
|
---|
1920 |
|
---|
1921 | VWN4LCFunctional::~VWN4LCFunctional()
|
---|
1922 | {
|
---|
1923 | }
|
---|
1924 |
|
---|
1925 | void
|
---|
1926 | VWN4LCFunctional::save_data_state(StateOut& s)
|
---|
1927 | {
|
---|
1928 | VWNLCFunctional::save_data_state(s);
|
---|
1929 | s.put(monte_carlo_prefactor_);
|
---|
1930 | }
|
---|
1931 |
|
---|
1932 | // based on the equations given on a NIST WWW site
|
---|
1933 | void
|
---|
1934 | VWN4LCFunctional::point_lc(const PointInputData &id, PointOutputData &od,
|
---|
1935 | double &ec_local, double &decrs, double &deczeta)
|
---|
1936 | {
|
---|
1937 | od.zero();
|
---|
1938 | const double fpp0 = 4./9. * 1./(pow(2., (1./3.)) - 1.);
|
---|
1939 | const double sixth = 1./6.;
|
---|
1940 | const double four_thirds = 4./3.;
|
---|
1941 | const double one_third = 1./3.;
|
---|
1942 | const double two_thirds = 2./3.;
|
---|
1943 |
|
---|
1944 | double rho = id.a.rho + id.b.rho;
|
---|
1945 | double zeta = (id.a.rho - id.b.rho)/rho;
|
---|
1946 | double x = pow(3./(4.*M_PI*rho), sixth);
|
---|
1947 | double rs = x*x;
|
---|
1948 |
|
---|
1949 | // Monte Carlo fitting parameters
|
---|
1950 | double epc_mc = F(x, Ap_, x0p_mc_, bp_mc_, cp_mc_);
|
---|
1951 | double efc_mc = F(x, Af_, x0f_mc_, bf_mc_, cf_mc_);
|
---|
1952 | // double alphac_mc = F(x, A_alpha_, x0_alpha_mc_, b_alpha_mc_, c_alpha_mc_);
|
---|
1953 | // RPA fitting parameters
|
---|
1954 | double epc_rpa = F(x, Ap_, x0p_rpa_, bp_rpa_, cp_rpa_);
|
---|
1955 | double efc_rpa = F(x, Af_, x0f_rpa_, bf_rpa_, cf_rpa_);
|
---|
1956 | double alphac_rpa = F(x, A_alpha_, x0_alpha_rpa_, b_alpha_rpa_, c_alpha_rpa_);
|
---|
1957 |
|
---|
1958 | double f = 9./8.*fpp0*(pow(1.+zeta, four_thirds)+pow(1.-zeta, four_thirds)-2.);
|
---|
1959 | double zeta2 = zeta*zeta;
|
---|
1960 | double zeta4 = zeta2*zeta2;
|
---|
1961 | double delta_e_rpa = efc_rpa - epc_rpa;
|
---|
1962 | double delta_e_mc = efc_mc - epc_mc;
|
---|
1963 | // double beta = fpp0 * delta_e_mc / alphac_rpa - 1.;
|
---|
1964 | // double delta_erpa_rszeta = alphac_rpa * f / fpp0 * (1. + beta * zeta4);
|
---|
1965 | // use delta_e_rpa here instead of delta_e_mc to get VWN3
|
---|
1966 | double delta_erpa_rszeta = alphac_rpa * f / fpp0 * (1. - zeta4)
|
---|
1967 | + f * delta_e_mc * zeta4;
|
---|
1968 | double delta_ec;
|
---|
1969 | if (!monte_carlo_prefactor_) delta_ec = delta_erpa_rszeta;
|
---|
1970 | else delta_ec = delta_e_mc/delta_e_rpa * delta_erpa_rszeta;
|
---|
1971 | // double ec = epc_rpa + delta_ec;
|
---|
1972 | double ec = epc_mc + delta_ec;
|
---|
1973 |
|
---|
1974 | od.energy = ec * rho;
|
---|
1975 | ec_local = ec;
|
---|
1976 |
|
---|
1977 | if (compute_potential_) {
|
---|
1978 | double zeta3 = zeta2*zeta;
|
---|
1979 | // Monte Carlo fitting parameters
|
---|
1980 | double depc_dr_s0_mc = dFdr_s(x, Ap_, x0p_mc_, bp_mc_, cp_mc_);
|
---|
1981 | double defc_dr_s1_mc = dFdr_s(x, Af_, x0f_mc_, bf_mc_, cf_mc_);
|
---|
1982 | // double dalphac_dr_s_mc = dFdr_s(x, A_alpha_, x0_alpha_mc_, b_alpha_mc_, c_alpha_mc_);
|
---|
1983 | // RPA fitting parameters
|
---|
1984 | double depc_dr_s0_rpa = dFdr_s(x, Ap_, x0p_rpa_, bp_rpa_, cp_rpa_);
|
---|
1985 | double defc_dr_s1_rpa = dFdr_s(x, Af_, x0f_rpa_, bf_rpa_, cf_rpa_);
|
---|
1986 | double dalphac_dr_s_rpa = dFdr_s(x, A_alpha_, x0_alpha_rpa_, b_alpha_rpa_,
|
---|
1987 | c_alpha_rpa_);
|
---|
1988 |
|
---|
1989 | double fp = two_thirds * (pow((1+zeta),one_third)
|
---|
1990 | - pow((1-zeta),one_third))/(pow(2.,one_third)-1);
|
---|
1991 |
|
---|
1992 | double ddelta_e_rpa = defc_dr_s1_rpa - depc_dr_s0_rpa;
|
---|
1993 | double ddelta_e_mc = defc_dr_s1_mc - depc_dr_s0_mc;
|
---|
1994 | // use ddelta_e_rpa here instead of ddelta_e_mc to get VWN3
|
---|
1995 | double ddeltae_rpa_dr_s = f / fpp0 * (1 - zeta4)* dalphac_dr_s_rpa
|
---|
1996 | + f * zeta4 * ddelta_e_mc;
|
---|
1997 | // use ddelta_e_rpa here instead of ddelta_e_mc to get VWN3
|
---|
1998 | double ddeltae_rpa_dzeta = alphac_rpa / fpp0 *
|
---|
1999 | ( fp * (1.-zeta4) - 4.* f * zeta3)
|
---|
2000 | + delta_e_mc * ( fp*zeta4 + 4.*f*zeta3);
|
---|
2001 |
|
---|
2002 | double dec_dzeta, dec_dr_s;
|
---|
2003 | if (!monte_carlo_prefactor_) {
|
---|
2004 | // dec_dr_s = depc_dr_s0_rpa + ddeltae_rpa_dr_s;
|
---|
2005 | dec_dr_s = depc_dr_s0_mc + ddeltae_rpa_dr_s;
|
---|
2006 | dec_dzeta = ddeltae_rpa_dzeta;
|
---|
2007 | }
|
---|
2008 | else {
|
---|
2009 | dec_dr_s = depc_dr_s0_rpa + delta_erpa_rszeta/delta_e_rpa * ddelta_e_mc
|
---|
2010 | + delta_e_mc/delta_e_rpa * ddeltae_rpa_dr_s
|
---|
2011 | - delta_erpa_rszeta*delta_e_mc/(delta_e_rpa*delta_e_rpa)* ddelta_e_rpa;
|
---|
2012 | dec_dzeta = delta_e_mc / delta_e_rpa * ddeltae_rpa_dzeta;
|
---|
2013 | }
|
---|
2014 |
|
---|
2015 | od.df_drho_a = ec - (rs/3.)*dec_dr_s - (zeta-1.)*dec_dzeta;
|
---|
2016 | od.df_drho_b = ec - (rs/3.)*dec_dr_s - (zeta+1.)*dec_dzeta;
|
---|
2017 | decrs = dec_dr_s;
|
---|
2018 | deczeta = dec_dzeta;
|
---|
2019 | }
|
---|
2020 | }
|
---|
2021 |
|
---|
2022 | /////////////////////////////////////////////////////////////////////////////
|
---|
2023 | // VWN5LCFunctional
|
---|
2024 | // Coded by Matt Leininger
|
---|
2025 |
|
---|
2026 | static ClassDesc VWN5LCFunctional_cd(
|
---|
2027 | typeid(VWN5LCFunctional),"VWN5LCFunctional",1,"public VWNLCFunctional",
|
---|
2028 | 0, create<VWN5LCFunctional>, create<VWN5LCFunctional>);
|
---|
2029 |
|
---|
2030 | VWN5LCFunctional::VWN5LCFunctional(StateIn& s):
|
---|
2031 | SavableState(s),
|
---|
2032 | VWNLCFunctional(s)
|
---|
2033 | {
|
---|
2034 | }
|
---|
2035 |
|
---|
2036 | VWN5LCFunctional::VWN5LCFunctional()
|
---|
2037 | {
|
---|
2038 | }
|
---|
2039 |
|
---|
2040 | VWN5LCFunctional::VWN5LCFunctional(const Ref<KeyVal>& keyval):
|
---|
2041 | VWNLCFunctional(keyval)
|
---|
2042 | {
|
---|
2043 | }
|
---|
2044 |
|
---|
2045 | VWN5LCFunctional::~VWN5LCFunctional()
|
---|
2046 | {
|
---|
2047 | }
|
---|
2048 |
|
---|
2049 | void
|
---|
2050 | VWN5LCFunctional::save_data_state(StateOut& s)
|
---|
2051 | {
|
---|
2052 | VWNLCFunctional::save_data_state(s);
|
---|
2053 | }
|
---|
2054 |
|
---|
2055 | // based on the equations given on a NIST WWW site
|
---|
2056 | // based on the VWN5 functional in Vosko, Wilk, and Nusair, Can. J. Phys.
|
---|
2057 | // 58, 1200, (1980). and Perdew and Wang, Phys. Rev. B, 45, 13244, (1992).
|
---|
2058 | void
|
---|
2059 | VWN5LCFunctional::point_lc(const PointInputData &id, PointOutputData &od,
|
---|
2060 | double &ec_local, double &decrs, double &deczeta)
|
---|
2061 | {
|
---|
2062 | od.zero();
|
---|
2063 | const double fpp0 = 4./9. * 1./(pow(2., (1./3.)) - 1.);
|
---|
2064 | const double sixth = 1./6.;
|
---|
2065 | const double four_thirds = 4./3.;
|
---|
2066 | const double one_third = 1./3.;
|
---|
2067 | const double two_thirds = 2./3.;
|
---|
2068 |
|
---|
2069 | double rho = id.a.rho + id.b.rho;
|
---|
2070 | double zeta = (id.a.rho - id.b.rho)/rho;
|
---|
2071 | double x = pow(3./(4.*M_PI*rho), sixth);
|
---|
2072 | double rs = x*x;
|
---|
2073 |
|
---|
2074 | double epc = F(x, Ap_, x0p_mc_, bp_mc_, cp_mc_);
|
---|
2075 | double efc = F(x, Af_, x0f_mc_, bf_mc_, cf_mc_);
|
---|
2076 | double alphac = F(x, A_alpha_, x0_alpha_mc_, b_alpha_mc_, c_alpha_mc_);
|
---|
2077 |
|
---|
2078 | double f = 9./8.*fpp0*(pow(1.+zeta, four_thirds)+pow(1.-zeta, four_thirds)-2.);
|
---|
2079 | double zeta2 = zeta*zeta;
|
---|
2080 | double zeta4 = zeta2*zeta2;
|
---|
2081 | double beta = fpp0 * (efc - epc) / alphac - 1.;
|
---|
2082 | double delta_ec = alphac * f / fpp0 * (1. + beta * zeta4);
|
---|
2083 | double ec = epc + delta_ec;
|
---|
2084 |
|
---|
2085 | ec_local = ec;
|
---|
2086 | od.energy = ec * rho;
|
---|
2087 |
|
---|
2088 | if (compute_potential_) {
|
---|
2089 | if (!spin_polarized_) {
|
---|
2090 | double depc_dr_s0 = dFdr_s(x, Ap_, x0p_mc_, bp_mc_, cp_mc_);
|
---|
2091 | double dec_dr_s = depc_dr_s0;
|
---|
2092 | od.df_drho_a = od.df_drho_b = ec - (rs/3.)*dec_dr_s;
|
---|
2093 | decrs = dec_dr_s;
|
---|
2094 | deczeta = 0.;
|
---|
2095 | }
|
---|
2096 | else {
|
---|
2097 | double zeta3 = zeta2*zeta;
|
---|
2098 | double depc_dr_s0 = dFdr_s(x, Ap_, x0p_mc_, bp_mc_, cp_mc_);
|
---|
2099 | double defc_dr_s1 = dFdr_s(x, Af_, x0f_mc_, bf_mc_, cf_mc_);
|
---|
2100 | double dalphac_dr_s = dFdr_s(x, A_alpha_, x0_alpha_mc_, b_alpha_mc_, c_alpha_mc_);
|
---|
2101 | double dec_dr_s = depc_dr_s0*(1 - f*zeta4) + defc_dr_s1 * f * zeta4
|
---|
2102 | + dalphac_dr_s * f / fpp0 * (1 - zeta4);
|
---|
2103 | double fp = two_thirds * (pow((1+zeta),one_third)
|
---|
2104 | - pow((1-zeta),one_third))/(pow(2.,one_third)-1);
|
---|
2105 | double dec_dzeta = 4.* zeta3 * f * (efc - epc - (alphac/fpp0))
|
---|
2106 | + fp * (zeta4 * (efc - epc) + (1-zeta4)*(alphac/fpp0));
|
---|
2107 | od.df_drho_a = ec - (rs/3.)*dec_dr_s - (zeta-1)*dec_dzeta;
|
---|
2108 | od.df_drho_b = ec - (rs/3.)*dec_dr_s - (zeta+1)*dec_dzeta;
|
---|
2109 | decrs = dec_dr_s;
|
---|
2110 | deczeta = dec_dzeta;
|
---|
2111 | }
|
---|
2112 | }
|
---|
2113 | }
|
---|
2114 |
|
---|
2115 | /////////////////////////////////////////////////////////////////////////////
|
---|
2116 | // XalphaFunctional
|
---|
2117 |
|
---|
2118 | static ClassDesc XalphaFunctional_cd(
|
---|
2119 | typeid(XalphaFunctional),"XalphaFunctional",1,"public DenFunctional",
|
---|
2120 | 0, create<XalphaFunctional>, create<XalphaFunctional>);
|
---|
2121 |
|
---|
2122 | XalphaFunctional::XalphaFunctional(StateIn& s):
|
---|
2123 | SavableState(s),
|
---|
2124 | DenFunctional(s)
|
---|
2125 | {
|
---|
2126 | s.get(alpha_);
|
---|
2127 | factor_ = alpha_ * 2.25 * pow(3.0/(4.*M_PI), 1.0/3.0);
|
---|
2128 | }
|
---|
2129 |
|
---|
2130 | XalphaFunctional::XalphaFunctional()
|
---|
2131 | {
|
---|
2132 | alpha_ = 0.70;
|
---|
2133 | factor_ = alpha_ * 2.25 * pow(3.0/(4.*M_PI), 1.0/3.0);
|
---|
2134 | }
|
---|
2135 |
|
---|
2136 | XalphaFunctional::XalphaFunctional(const Ref<KeyVal>& keyval):
|
---|
2137 | DenFunctional(keyval)
|
---|
2138 | {
|
---|
2139 | alpha_ = keyval->doublevalue("alpha", KeyValValuedouble(0.70));
|
---|
2140 | factor_ = alpha_ * 2.25 * pow(3.0/(4.*M_PI), 1.0/3.0);
|
---|
2141 | }
|
---|
2142 |
|
---|
2143 | XalphaFunctional::~XalphaFunctional()
|
---|
2144 | {
|
---|
2145 | }
|
---|
2146 |
|
---|
2147 | void
|
---|
2148 | XalphaFunctional::save_data_state(StateOut& s)
|
---|
2149 | {
|
---|
2150 | DenFunctional::save_data_state(s);
|
---|
2151 | s.put(alpha_);
|
---|
2152 | }
|
---|
2153 |
|
---|
2154 | void
|
---|
2155 | XalphaFunctional::point(const PointInputData &id,
|
---|
2156 | PointOutputData &od)
|
---|
2157 | {
|
---|
2158 | od.zero();
|
---|
2159 | const double four_thirds = 4./3.;
|
---|
2160 |
|
---|
2161 | if (!spin_polarized_) {
|
---|
2162 | od.energy = - 2.0 * factor_ * id.a.rho * id.a.rho_13;
|
---|
2163 | if (compute_potential_) {
|
---|
2164 | od.df_drho_a = -four_thirds * factor_ * id.a.rho_13;
|
---|
2165 | od.df_drho_b = od.df_drho_a;
|
---|
2166 | }
|
---|
2167 | }
|
---|
2168 | else {
|
---|
2169 | double rhoa43 = id.a.rho * id.a.rho_13;
|
---|
2170 | double rhob43 = id.b.rho * id.b.rho_13;
|
---|
2171 | od.energy = - factor_ * (rhoa43 + rhob43);
|
---|
2172 | if (compute_potential_) {
|
---|
2173 | od.df_drho_a = -four_thirds * factor_ * id.a.rho_13;
|
---|
2174 | od.df_drho_b = -four_thirds * factor_ * id.b.rho_13;
|
---|
2175 | }
|
---|
2176 | }
|
---|
2177 | }
|
---|
2178 |
|
---|
2179 | void
|
---|
2180 | XalphaFunctional::print(ostream& o) const
|
---|
2181 | {
|
---|
2182 | o
|
---|
2183 | << indent << scprintf("XalphaFunctional: alpha = %12.8f", alpha_) << endl;
|
---|
2184 | }
|
---|
2185 |
|
---|
2186 | /////////////////////////////////////////////////////////////////////////////
|
---|
2187 | // Becke88XFunctional
|
---|
2188 |
|
---|
2189 | static ClassDesc Becke88XFunctional_cd(
|
---|
2190 | typeid(Becke88XFunctional),"Becke88XFunctional",1,"public DenFunctional",
|
---|
2191 | 0, create<Becke88XFunctional>, create<Becke88XFunctional>);
|
---|
2192 |
|
---|
2193 | Becke88XFunctional::Becke88XFunctional(StateIn& s):
|
---|
2194 | SavableState(s),
|
---|
2195 | DenFunctional(s)
|
---|
2196 | {
|
---|
2197 | s.get(beta_);
|
---|
2198 | beta6_ = 6. * beta_;
|
---|
2199 | }
|
---|
2200 |
|
---|
2201 | Becke88XFunctional::Becke88XFunctional()
|
---|
2202 | {
|
---|
2203 | beta_ = 0.0042;
|
---|
2204 | beta6_ = 6. * beta_;
|
---|
2205 | }
|
---|
2206 |
|
---|
2207 | Becke88XFunctional::Becke88XFunctional(const Ref<KeyVal>& keyval):
|
---|
2208 | DenFunctional(keyval)
|
---|
2209 | {
|
---|
2210 | beta_ = keyval->doublevalue("beta", KeyValValuedouble(0.0042));
|
---|
2211 | beta6_ = 6. * beta_;
|
---|
2212 | }
|
---|
2213 |
|
---|
2214 | Becke88XFunctional::~Becke88XFunctional()
|
---|
2215 | {
|
---|
2216 | }
|
---|
2217 |
|
---|
2218 | void
|
---|
2219 | Becke88XFunctional::save_data_state(StateOut& s)
|
---|
2220 | {
|
---|
2221 | DenFunctional::save_data_state(s);
|
---|
2222 | s.put(beta_);
|
---|
2223 | }
|
---|
2224 |
|
---|
2225 | int
|
---|
2226 | Becke88XFunctional::need_density_gradient()
|
---|
2227 | {
|
---|
2228 | return 1;
|
---|
2229 | }
|
---|
2230 |
|
---|
2231 | // Becke's exchange
|
---|
2232 | // From: C.W. Murray et al. Mol. Phys. Vol 78 pp 997-1014 (1993)
|
---|
2233 | // originally coded by Mike Colvin
|
---|
2234 | void
|
---|
2235 | Becke88XFunctional::point(const PointInputData &id,
|
---|
2236 | PointOutputData &od)
|
---|
2237 | {
|
---|
2238 | od.zero();
|
---|
2239 | // Preset terms from Murray's paper
|
---|
2240 | const double beta = beta_;
|
---|
2241 | const double beta6 = beta6_;
|
---|
2242 |
|
---|
2243 | // const double beta6=0.0252;
|
---|
2244 |
|
---|
2245 | // Use simplified formula
|
---|
2246 | double rho_a_13 = pow(id.a.rho,(1./3.));
|
---|
2247 | double rho_a_43 = id.a.rho*rho_a_13;
|
---|
2248 | double xa, xa2, ga_denom, ga_denom2;
|
---|
2249 | if (id.a.rho > MIN_DENSITY) {
|
---|
2250 | xa = sqrt(id.a.gamma)/rho_a_43;
|
---|
2251 | xa2 = xa*xa;
|
---|
2252 | ga_denom = 1/(1.+beta6*xa*asinh(xa));
|
---|
2253 | ga_denom2 = ga_denom*ga_denom;
|
---|
2254 | }
|
---|
2255 | else {
|
---|
2256 | xa = xa2 = 0.;
|
---|
2257 | ga_denom = ga_denom2 = 1.;
|
---|
2258 | }
|
---|
2259 | double Fa = sqrt(1.+xa2);
|
---|
2260 | double Ha = 1. - 6.*beta*xa2/Fa;
|
---|
2261 | double ex;
|
---|
2262 | if (id.a.rho > MIN_DENSITY) ex = -rho_a_43*beta*xa2*ga_denom;
|
---|
2263 | else ex = 0.;
|
---|
2264 |
|
---|
2265 | if (compute_potential_) {
|
---|
2266 | if (id.a.rho > MIN_DENSITY) {
|
---|
2267 | od.df_drho_a = 4./3. * beta * rho_a_13 * xa2 * ga_denom2 * Ha;
|
---|
2268 | od.df_dgamma_aa = -beta * ga_denom / (2.*rho_a_43) * (1. + ga_denom*Ha);
|
---|
2269 | }
|
---|
2270 | else od.df_drho_a = od.df_dgamma_aa = 0.;
|
---|
2271 | od.df_drho_b=od.df_drho_a;
|
---|
2272 | od.df_dgamma_bb=od.df_dgamma_aa;
|
---|
2273 | }
|
---|
2274 |
|
---|
2275 | if (spin_polarized_) {
|
---|
2276 | double rho_b_13 = pow(id.b.rho,(1./3.));
|
---|
2277 | double rho_b_43 = id.b.rho*rho_b_13;
|
---|
2278 | double xb, xb2, gb_denom, gb_denom2;
|
---|
2279 | if (id.b.rho > MIN_DENSITY) {
|
---|
2280 | xb = sqrt(id.b.gamma)/rho_b_43;
|
---|
2281 | xb2 = xb*xb;
|
---|
2282 | gb_denom = 1./(1.+beta6*xb*asinh(xb));
|
---|
2283 | gb_denom2 = gb_denom*gb_denom;
|
---|
2284 | }
|
---|
2285 | else {
|
---|
2286 | xb = xb2 = 0.;
|
---|
2287 | gb_denom = gb_denom2 = 1.;
|
---|
2288 | }
|
---|
2289 | double Fb = sqrt(1.+xb2);
|
---|
2290 | double Hb = 1. - 6.*beta*xb2/Fb;
|
---|
2291 | ex += -rho_b_43*beta*xb2*gb_denom;
|
---|
2292 |
|
---|
2293 | if (compute_potential_) {
|
---|
2294 | if (id.b.rho > MIN_DENSITY) {
|
---|
2295 | od.df_drho_b = 4./3. * beta * rho_b_13 * xb2 * gb_denom2 * Hb;
|
---|
2296 | od.df_dgamma_bb = -beta / (2.*rho_b_43) * (gb_denom + gb_denom2*Hb);
|
---|
2297 | }
|
---|
2298 | else od.df_drho_b = od.df_dgamma_bb = 0.;
|
---|
2299 | }
|
---|
2300 | }
|
---|
2301 | else ex += ex;
|
---|
2302 |
|
---|
2303 | od.energy = ex;
|
---|
2304 | }
|
---|
2305 |
|
---|
2306 |
|
---|
2307 | /////////////////////////////////////////////////////////////////////////////
|
---|
2308 | // LYPCFunctional
|
---|
2309 | // Coded by Matt Leininger
|
---|
2310 |
|
---|
2311 | static ClassDesc LYPCFunctional_cd(
|
---|
2312 | typeid(LYPCFunctional),"LYPCFunctional",1,"public DenFunctional",
|
---|
2313 | 0, create<LYPCFunctional>, create<LYPCFunctional>);
|
---|
2314 |
|
---|
2315 | LYPCFunctional::LYPCFunctional(StateIn& s):
|
---|
2316 | SavableState(s),
|
---|
2317 | DenFunctional(s)
|
---|
2318 | {
|
---|
2319 | s.get(a_);
|
---|
2320 | s.get(b_);
|
---|
2321 | s.get(c_);
|
---|
2322 | s.get(d_);
|
---|
2323 | }
|
---|
2324 |
|
---|
2325 | LYPCFunctional::LYPCFunctional()
|
---|
2326 | {
|
---|
2327 | init_constants();
|
---|
2328 | }
|
---|
2329 |
|
---|
2330 | LYPCFunctional::LYPCFunctional(const Ref<KeyVal>& keyval):
|
---|
2331 | DenFunctional(keyval)
|
---|
2332 | {
|
---|
2333 | init_constants();
|
---|
2334 | a_ = keyval->doublevalue("a", KeyValValuedouble(a_));
|
---|
2335 | b_ = keyval->doublevalue("b", KeyValValuedouble(b_));
|
---|
2336 | c_ = keyval->doublevalue("c", KeyValValuedouble(c_));
|
---|
2337 | d_ = keyval->doublevalue("d", KeyValValuedouble(d_));
|
---|
2338 | }
|
---|
2339 |
|
---|
2340 | LYPCFunctional::~LYPCFunctional()
|
---|
2341 | {
|
---|
2342 | }
|
---|
2343 |
|
---|
2344 | void
|
---|
2345 | LYPCFunctional::save_data_state(StateOut& s)
|
---|
2346 | {
|
---|
2347 | DenFunctional::save_data_state(s);
|
---|
2348 | s.put(a_);
|
---|
2349 | s.put(b_);
|
---|
2350 | s.put(c_);
|
---|
2351 | s.put(d_);
|
---|
2352 | }
|
---|
2353 |
|
---|
2354 | void
|
---|
2355 | LYPCFunctional::init_constants()
|
---|
2356 | {
|
---|
2357 | a_ = 0.04918;
|
---|
2358 | b_ = 0.132;
|
---|
2359 | c_ = 0.2533;
|
---|
2360 | d_ = 0.349;
|
---|
2361 | }
|
---|
2362 |
|
---|
2363 | int
|
---|
2364 | LYPCFunctional::need_density_gradient()
|
---|
2365 | {
|
---|
2366 | return 1;
|
---|
2367 | }
|
---|
2368 |
|
---|
2369 | // Lee-Yang-Parr correlation
|
---|
2370 | // From: Burkhard Miehlich, et al. Chem Phys. Lett. Vol. 157 pp200-206 (1989)
|
---|
2371 | // Original LYP paper Phys. Rev. B Vol. 37 pp785-789 (1988)
|
---|
2372 | // originally coded by Mike Colvin
|
---|
2373 | // potential terms added by Matt Leininger
|
---|
2374 | void
|
---|
2375 | LYPCFunctional::point(const PointInputData &id,
|
---|
2376 | PointOutputData &od)
|
---|
2377 | {
|
---|
2378 | od.zero();
|
---|
2379 | double ec;
|
---|
2380 |
|
---|
2381 | // Precalculate terms for efficiency
|
---|
2382 | double dens=id.a.rho+id.b.rho;
|
---|
2383 | double dens2=dens*dens;
|
---|
2384 | // double grad=id.a.gamma+id.b.gamma;
|
---|
2385 | // double grad=sqrt( pow(id.a.del_rho[0]+id.b.del_rho[0],2.)
|
---|
2386 | // + pow(id.a.del_rho[1]+id.b.del_rho[1],2.)
|
---|
2387 | // + pow(id.a.del_rho[2]+id.b.del_rho[2],2.) );
|
---|
2388 | // double grad2=grad*grad;
|
---|
2389 | double dens1_3=pow(dens,-1./3.);
|
---|
2390 |
|
---|
2391 | // Precalculate terms defined in Miehlich's paper
|
---|
2392 | const double a = a_;
|
---|
2393 | const double b = b_;
|
---|
2394 | const double c = c_;
|
---|
2395 | const double d = d_;
|
---|
2396 | double omega=exp(-c*dens1_3)/(1.+d*dens1_3)*pow(dens,-11./3.);
|
---|
2397 | double delta=c*dens1_3+d*dens1_3/(1.+d*dens1_3);
|
---|
2398 | double cf=0.3*pow(3.* M_PI*M_PI,2./3.);
|
---|
2399 | double denom = 1.+d*dens1_3;
|
---|
2400 |
|
---|
2401 | double dens_a2=id.a.rho*id.a.rho;
|
---|
2402 | double dens_b2=id.b.rho*id.b.rho;
|
---|
2403 | double dens_ab=id.a.rho*id.b.rho;
|
---|
2404 | double grad_a2=id.a.gamma;
|
---|
2405 | double grad_b2=id.b.gamma;
|
---|
2406 | double grad_ab=id.gamma_ab;
|
---|
2407 |
|
---|
2408 | double eflyp_1 = -4.*a*dens_ab/(dens*denom);
|
---|
2409 | double intermediate_1 = pow(2.,2./3.)*144.*cf*(pow(id.a.rho,8./3.)+pow(id.b.rho,8./3.))
|
---|
2410 | + (47.-7.*delta)*(grad_a2+grad_b2+2.*grad_ab)-(45.-delta)*(grad_a2+grad_b2)
|
---|
2411 | + 2.*(11.-delta)/dens*(id.a.rho*grad_a2+id.b.rho*grad_b2);
|
---|
2412 | double intermediate_2 = -4./3.*dens2*grad_ab - (dens_a2*grad_b2+dens_b2*grad_a2);
|
---|
2413 | double intermediate_3 = dens_ab/18.* intermediate_1 + intermediate_2;
|
---|
2414 | double eflyp_2 = -omega*a*b*intermediate_3;
|
---|
2415 | ec = eflyp_1 + eflyp_2;
|
---|
2416 | od.energy = ec;
|
---|
2417 |
|
---|
2418 | if (compute_potential_) {
|
---|
2419 | double dens4_3 = pow(dens,-4./3.);
|
---|
2420 | double expo = exp(-c*dens1_3);
|
---|
2421 | double ddelta_drho_a = 1./3* (d*d*dens4_3*dens1_3/(denom*denom)
|
---|
2422 | - delta/dens);
|
---|
2423 | double domega_drho_a = -1./3.*omega*dens4_3*(11./dens1_3 - c - d/denom);
|
---|
2424 | double df1_drho_a;
|
---|
2425 | df1_drho_a = -4.*a*id.b.rho/(dens*denom) *
|
---|
2426 | (id.a.rho/3.*d*dens4_3/denom + 1. - id.a.rho/dens);
|
---|
2427 | // if (id.a.rho > MIN_DENSITY)
|
---|
2428 | // df1_drho_a = -4.*a*dens_ab/(dens*denom)
|
---|
2429 | // * (1./3.*d*dens4_3/denom + 1/id.a.rho - 1./dens);
|
---|
2430 | // else df1_drho_a = 0.;
|
---|
2431 | double df2_drho_a = -domega_drho_a*a*b*intermediate_3
|
---|
2432 | - omega*a*b*( id.b.rho/18.* intermediate_1
|
---|
2433 | + dens_ab/18.*(144.*pow(2.,2./3.)*cf*8./3.*pow(id.a.rho,5./3.)
|
---|
2434 | + 2.*(11.-delta)*grad_a2/dens
|
---|
2435 | - 2./dens*ddelta_drho_a*(id.a.rho*grad_a2 + id.b.rho*grad_b2)
|
---|
2436 | - 2.*(11.-delta)/dens2*(id.a.rho*grad_a2 + id.b.rho*grad_b2)
|
---|
2437 | - 7.*ddelta_drho_a*(grad_a2+grad_b2+2.*grad_ab)
|
---|
2438 | + ddelta_drho_a*(grad_a2+grad_b2) )
|
---|
2439 | - 8./3.*dens*grad_ab - 2.*id.a.rho*grad_b2 );
|
---|
2440 |
|
---|
2441 | od.df_drho_a = df1_drho_a + df2_drho_a;
|
---|
2442 |
|
---|
2443 | od.df_dgamma_aa = -omega*a*b
|
---|
2444 | * (dens_ab/9.*(1.-3.*delta + id.a.rho*(11.-delta)/dens) - dens_b2);
|
---|
2445 | od.df_dgamma_ab = -omega*a*b*(dens_ab/9.*(47.-7.*delta) - 4./3.*dens2);
|
---|
2446 | od.df_drho_b = od.df_drho_a;
|
---|
2447 | od.df_dgamma_bb = od.df_dgamma_aa;
|
---|
2448 | if (spin_polarized_) {
|
---|
2449 | double ddelta_drho_b = ddelta_drho_a;
|
---|
2450 | double domega_drho_b = domega_drho_a;
|
---|
2451 | double df1_drho_b;
|
---|
2452 | df1_drho_b = -4.*a*id.a.rho/(dens*denom) *
|
---|
2453 | (id.b.rho/3.*d*dens4_3/denom + 1. - id.b.rho/dens);
|
---|
2454 | //if (id.b.rho > MIN_DENSITY)
|
---|
2455 | // df1_drho_b = -4.*a*dens_ab/(dens*denom)
|
---|
2456 | // * (1./3.*d*dens4_3/denom + 1./id.b.rho - 1./dens);
|
---|
2457 | //else df1_drho_b = 0.;
|
---|
2458 | double df2_drho_b = -domega_drho_b*a*b*intermediate_3
|
---|
2459 | - omega*a*b*( id.a.rho/18.* intermediate_1
|
---|
2460 | + dens_ab/18.*(144.*pow(2.,2./3.)*cf*8./3.*pow(id.b.rho,5./3.)
|
---|
2461 | + 2.*(11.-delta)*grad_b2/dens
|
---|
2462 | - 2./dens*ddelta_drho_b*(id.a.rho*grad_a2 + id.b.rho*grad_b2)
|
---|
2463 | - 2.*(11.-delta)/dens2*(id.a.rho*grad_a2 + id.b.rho*grad_b2)
|
---|
2464 | - 7.*ddelta_drho_b*(grad_a2+grad_b2+2.*grad_ab)
|
---|
2465 | + ddelta_drho_b*(grad_a2+grad_b2) )
|
---|
2466 | - 8./3.*dens*grad_ab - 2.*id.b.rho*grad_a2 );
|
---|
2467 | od.df_drho_b = df1_drho_b + df2_drho_b;
|
---|
2468 | od.df_dgamma_bb = -omega*a*b
|
---|
2469 | * (dens_ab/9.*(1.-3.*delta + id.b.rho*(11.-delta)/dens) - dens_a2);
|
---|
2470 | }
|
---|
2471 | }
|
---|
2472 |
|
---|
2473 | }
|
---|
2474 |
|
---|
2475 | /////////////////////////////////////////////////////////////////////////////
|
---|
2476 | // Perdew 1986 (P86) Correlation Functional
|
---|
2477 | // J. P. Perdew, PRB, 33, 8822, 1986.
|
---|
2478 | // C. W. Murray, N. C. Handy and G. J. Laming, Mol. Phys., 78, 997, 1993.
|
---|
2479 | //
|
---|
2480 | // Coded by Matt Leininger
|
---|
2481 |
|
---|
2482 | static ClassDesc P86CFunctional_cd(
|
---|
2483 | typeid(P86CFunctional),"P86CFunctional",1,"public DenFunctional",
|
---|
2484 | 0, create<P86CFunctional>, create<P86CFunctional>);
|
---|
2485 |
|
---|
2486 | P86CFunctional::P86CFunctional(StateIn& s):
|
---|
2487 | SavableState(s),
|
---|
2488 | DenFunctional(s)
|
---|
2489 | {
|
---|
2490 | s.get(a_);
|
---|
2491 | s.get(C1_);
|
---|
2492 | s.get(C2_);
|
---|
2493 | s.get(C3_);
|
---|
2494 | s.get(C4_);
|
---|
2495 | s.get(C5_);
|
---|
2496 | s.get(C6_);
|
---|
2497 | s.get(C7_);
|
---|
2498 | }
|
---|
2499 |
|
---|
2500 | P86CFunctional::P86CFunctional()
|
---|
2501 | {
|
---|
2502 | init_constants();
|
---|
2503 | }
|
---|
2504 |
|
---|
2505 | P86CFunctional::P86CFunctional(const Ref<KeyVal>& keyval):
|
---|
2506 | DenFunctional(keyval)
|
---|
2507 | {
|
---|
2508 | init_constants();
|
---|
2509 | a_ = keyval->doublevalue("a", KeyValValuedouble(a_));
|
---|
2510 | C1_ = keyval->doublevalue("C1", KeyValValuedouble(C1_));
|
---|
2511 | C2_ = keyval->doublevalue("C2", KeyValValuedouble(C2_));
|
---|
2512 | C3_ = keyval->doublevalue("C3", KeyValValuedouble(C3_));
|
---|
2513 | C4_ = keyval->doublevalue("C4", KeyValValuedouble(C4_));
|
---|
2514 | C5_ = keyval->doublevalue("C5", KeyValValuedouble(C5_));
|
---|
2515 | C6_ = keyval->doublevalue("C6", KeyValValuedouble(C6_));
|
---|
2516 | C7_ = keyval->doublevalue("C7", KeyValValuedouble(C7_));
|
---|
2517 | }
|
---|
2518 |
|
---|
2519 | P86CFunctional::~P86CFunctional()
|
---|
2520 | {
|
---|
2521 | }
|
---|
2522 |
|
---|
2523 | void
|
---|
2524 | P86CFunctional::save_data_state(StateOut& s)
|
---|
2525 | {
|
---|
2526 | DenFunctional::save_data_state(s);
|
---|
2527 | s.put(a_);
|
---|
2528 | s.put(C1_);
|
---|
2529 | s.put(C2_);
|
---|
2530 | s.put(C3_);
|
---|
2531 | s.put(C4_);
|
---|
2532 | s.put(C5_);
|
---|
2533 | s.put(C6_);
|
---|
2534 | s.put(C7_);
|
---|
2535 | }
|
---|
2536 |
|
---|
2537 | void
|
---|
2538 | P86CFunctional::init_constants()
|
---|
2539 | {
|
---|
2540 | a_ = 1.745*0.11;
|
---|
2541 | C1_ = 0.001667;
|
---|
2542 | C2_ = 0.002568;
|
---|
2543 | C3_ = 0.023266;
|
---|
2544 | C4_ = 7.389e-6;
|
---|
2545 | C5_ = 8.723;
|
---|
2546 | C6_ = 0.472;
|
---|
2547 | C7_ = 1e4*C4_;
|
---|
2548 | }
|
---|
2549 |
|
---|
2550 | int
|
---|
2551 | P86CFunctional::need_density_gradient()
|
---|
2552 | {
|
---|
2553 | return 1;
|
---|
2554 | }
|
---|
2555 |
|
---|
2556 | void
|
---|
2557 | P86CFunctional::point(const PointInputData &id,
|
---|
2558 | PointOutputData &od)
|
---|
2559 | {
|
---|
2560 | od.zero();
|
---|
2561 | // Precalculate terms for efficiency
|
---|
2562 | double rho = id.a.rho + id.b.rho;
|
---|
2563 | double rho76 = pow(rho, (7./6.));
|
---|
2564 | double rho43 = pow(rho, (4./3.));
|
---|
2565 | double rho13 = pow(rho, (1./3.));
|
---|
2566 | double rho16 = pow(rho, (1./6.));
|
---|
2567 | double rs = pow( (3./(4.*M_PI*rho)), (1./3.));
|
---|
2568 | double rs2 = rs*rs;
|
---|
2569 | double rs3 = rs2*rs;
|
---|
2570 | double zeta = (id.a.rho - id.b.rho)/rho;
|
---|
2571 | double fzeta = pow(2.,(1./3.)) *
|
---|
2572 | sqrt( pow(((1.+zeta)/2.),(5./3.)) + pow(((1.-zeta)/2.),(5./3.)) );
|
---|
2573 | double C_infin = C1_ + C2_;
|
---|
2574 | double numer = C2_ + C3_*rs + C4_*rs2;
|
---|
2575 | double denom = 1.+C5_*rs+C6_*rs2+C7_*rs3;
|
---|
2576 | double C_rho = C1_ + numer/denom;
|
---|
2577 | double gamma_aa = id.a.gamma;
|
---|
2578 | double gamma_bb = id.b.gamma;
|
---|
2579 | double gamma_ab = id.gamma_ab;
|
---|
2580 | double gamma_total = sqrt(gamma_aa + gamma_bb + 2.*gamma_ab);
|
---|
2581 | double gamma_total2 = gamma_total*gamma_total;
|
---|
2582 | double Phi = a_ * C_infin * gamma_total / (C_rho * rho76);
|
---|
2583 | double fp86 = exp(-Phi) * C_rho * gamma_total2 / (fzeta*rho43);
|
---|
2584 | if (rho < MIN_DENSITY) fp86 = 0.;
|
---|
2585 | od.energy = fp86;
|
---|
2586 |
|
---|
2587 | if (compute_potential_) {
|
---|
2588 | double drs_drhoa = -rs/(3.*rho);
|
---|
2589 | double dCrho_drhoa = drs_drhoa/denom *
|
---|
2590 | (C3_+2.*C4_*rs - numer/denom * (C5_+2.*C6_*rs+3.*C7_*rs2));
|
---|
2591 | double dCrho_drhob = dCrho_drhoa;
|
---|
2592 | double dzeta_drhoa = 1./rho * (1.-zeta);
|
---|
2593 | double dzeta_drhob = 1./rho * (-1.-zeta);
|
---|
2594 | double dPhi_drhoa = -Phi/(C_rho*rho76)*(dCrho_drhoa*rho76 + C_rho*(7./6.)*rho16);
|
---|
2595 | double dPhi_drhob = dPhi_drhoa;
|
---|
2596 | double dfzeta_drhoa = pow(2., (-1./3.))*1./fzeta *
|
---|
2597 | (5./3. * pow(((1.+zeta)/2.), (2./3.))*0.5*dzeta_drhoa
|
---|
2598 | + 5./3.*pow(((1.-zeta)/2.), (2./3.))*-0.5*dzeta_drhoa);
|
---|
2599 | double dfzeta_drhob = pow(2., (-1./3.))*1./fzeta *
|
---|
2600 | (5./3. * pow(((1.+zeta)/2.), (2./3.))*0.5*dzeta_drhob
|
---|
2601 | + 5./3.*pow(((1.-zeta)/2.), (2./3.))*-0.5*dzeta_drhob);
|
---|
2602 | double dfp86_drhoa = fp86/C_rho*(-dPhi_drhoa*C_rho + dCrho_drhoa)
|
---|
2603 | - fp86/(fzeta*rho43)*(dfzeta_drhoa*rho43 + fzeta*(4./3.)*rho13);
|
---|
2604 | double dfp86_drhob = fp86/C_rho*(-dPhi_drhob*C_rho + dCrho_drhob)
|
---|
2605 | - fp86/(fzeta*rho43)*(dfzeta_drhob*rho43 + fzeta*(4./3.)*rho13);
|
---|
2606 | if (rho < MIN_DENSITY) dfp86_drhoa = dfp86_drhob = 0.;
|
---|
2607 | od.df_drho_a = dfp86_drhoa;
|
---|
2608 | od.df_drho_b = dfp86_drhob;
|
---|
2609 |
|
---|
2610 | // gamma part of potential
|
---|
2611 | // double dPhi_dgamma_aa = Phi/(2.*gamma_total2);
|
---|
2612 | // double dPhi_dgamma_bb = dPhi_dgamma_aa;
|
---|
2613 | // double dfp86_dgamma_aa = fp86*(1./gamma_total2 - dPhi_dgamma_aa);
|
---|
2614 | double prefactor = exp(-Phi)*C_rho/( fzeta*rho43 );
|
---|
2615 | double dfp86_dgamma_aa = prefactor * (1.-Phi/2.);
|
---|
2616 | double dfp86_dgamma_bb = dfp86_dgamma_aa;
|
---|
2617 | if (rho < MIN_DENSITY) dfp86_dgamma_aa = dfp86_dgamma_bb = 0.;
|
---|
2618 | od.df_dgamma_aa = dfp86_dgamma_aa;
|
---|
2619 | od.df_dgamma_bb = dfp86_dgamma_bb;
|
---|
2620 |
|
---|
2621 | // double dPhi_dgamma_ab = 2.*dPhi_dgamma_aa;
|
---|
2622 | // double dfp86_dgamma_ab = fp86*(2./gamma_total2 - dPhi_dgamma_ab);
|
---|
2623 | double dfp86_dgamma_ab = prefactor * (2.-Phi);
|
---|
2624 | if (rho < MIN_DENSITY) dfp86_dgamma_ab= 0.;
|
---|
2625 | od.df_dgamma_ab = dfp86_dgamma_ab;
|
---|
2626 |
|
---|
2627 | }
|
---|
2628 | }
|
---|
2629 |
|
---|
2630 |
|
---|
2631 | /////////////////////////////////////////////////////////////////////////////
|
---|
2632 | // Perdew 1986 (P86) Correlation Functional
|
---|
2633 | // J. P. Perdew, PRB, 33, 8822, 1986.
|
---|
2634 | // C. W. Murray, N. C. Handy and G. J. Laming, Mol. Phys., 78, 997, 1993.
|
---|
2635 | //
|
---|
2636 | // Coded by Matt Leininger
|
---|
2637 |
|
---|
2638 | static ClassDesc NewP86CFunctional_cd(
|
---|
2639 | typeid(NewP86CFunctional),"NewP86CFunctional",1,"public DenFunctional",
|
---|
2640 | 0, create<NewP86CFunctional>, create<NewP86CFunctional>);
|
---|
2641 |
|
---|
2642 | NewP86CFunctional::NewP86CFunctional(StateIn& s):
|
---|
2643 | SavableState(s),
|
---|
2644 | DenFunctional(s)
|
---|
2645 | {
|
---|
2646 | s.get(a_);
|
---|
2647 | s.get(C1_);
|
---|
2648 | s.get(C2_);
|
---|
2649 | s.get(C3_);
|
---|
2650 | s.get(C4_);
|
---|
2651 | s.get(C5_);
|
---|
2652 | s.get(C6_);
|
---|
2653 | s.get(C7_);
|
---|
2654 | }
|
---|
2655 |
|
---|
2656 | NewP86CFunctional::NewP86CFunctional()
|
---|
2657 | {
|
---|
2658 | init_constants();
|
---|
2659 | }
|
---|
2660 |
|
---|
2661 | NewP86CFunctional::NewP86CFunctional(const Ref<KeyVal>& keyval):
|
---|
2662 | DenFunctional(keyval)
|
---|
2663 | {
|
---|
2664 | init_constants();
|
---|
2665 | a_ = keyval->doublevalue("a", KeyValValuedouble(a_));
|
---|
2666 | C1_ = keyval->doublevalue("C1", KeyValValuedouble(C1_));
|
---|
2667 | C2_ = keyval->doublevalue("C2", KeyValValuedouble(C2_));
|
---|
2668 | C3_ = keyval->doublevalue("C3", KeyValValuedouble(C3_));
|
---|
2669 | C4_ = keyval->doublevalue("C4", KeyValValuedouble(C4_));
|
---|
2670 | C5_ = keyval->doublevalue("C5", KeyValValuedouble(C5_));
|
---|
2671 | C6_ = keyval->doublevalue("C6", KeyValValuedouble(C6_));
|
---|
2672 | C7_ = keyval->doublevalue("C7", KeyValValuedouble(C7_));
|
---|
2673 | }
|
---|
2674 |
|
---|
2675 | NewP86CFunctional::~NewP86CFunctional()
|
---|
2676 | {
|
---|
2677 | }
|
---|
2678 |
|
---|
2679 | void
|
---|
2680 | NewP86CFunctional::save_data_state(StateOut& s)
|
---|
2681 | {
|
---|
2682 | DenFunctional::save_data_state(s);
|
---|
2683 | s.put(a_);
|
---|
2684 | s.put(C1_);
|
---|
2685 | s.put(C2_);
|
---|
2686 | s.put(C3_);
|
---|
2687 | s.put(C4_);
|
---|
2688 | s.put(C5_);
|
---|
2689 | s.put(C6_);
|
---|
2690 | s.put(C7_);
|
---|
2691 | }
|
---|
2692 |
|
---|
2693 | void
|
---|
2694 | NewP86CFunctional::init_constants()
|
---|
2695 | {
|
---|
2696 | a_ = 1.745*0.11;
|
---|
2697 | C1_ = 0.001667;
|
---|
2698 | C2_ = 0.002568;
|
---|
2699 | C3_ = 0.023266;
|
---|
2700 | C4_ = 7.389e-6;
|
---|
2701 | C5_ = 8.723;
|
---|
2702 | C6_ = 0.472;
|
---|
2703 | C7_ = 1e4*C4_;
|
---|
2704 | }
|
---|
2705 |
|
---|
2706 | int
|
---|
2707 | NewP86CFunctional::need_density_gradient()
|
---|
2708 | {
|
---|
2709 | return 1;
|
---|
2710 | }
|
---|
2711 |
|
---|
2712 | double
|
---|
2713 | NewP86CFunctional::rho_deriv(double rho_a, double rho_b, double mdr)
|
---|
2714 | {
|
---|
2715 | double ra0,ra1,ra2,ra3,ra4,ra5,ra6,ra7,ra8,ra9,ra10,ra11,ra12,ra13,ra14,
|
---|
2716 | ra15,ra16,ra17,ra18,ra19,ra20;
|
---|
2717 | double c_infin, c_1, c_2, c_3, c_4, c_5, c_6, c_7;
|
---|
2718 |
|
---|
2719 | c_infin = C1_ + C2_;
|
---|
2720 | c_1 = C1_; c_2 = C2_; c_3 = C3_; c_4 = C4_; c_5 = C5_; c_6 = C6_; c_7 = C7_;
|
---|
2721 |
|
---|
2722 | ra0 = pow(mdr,2.0);
|
---|
2723 | ra1 = rho_b+rho_a;
|
---|
2724 | ra2 = 1/pow(ra1,1.3333333333333333);
|
---|
2725 | ra3 = -(1.0*rho_b)+rho_a;
|
---|
2726 | ra4 = 1/ra1;
|
---|
2727 | // ra4 = 1/pow(ra1,1.0);
|
---|
2728 | ra5 = -(1.0*ra3*ra4)+1.0;
|
---|
2729 | ra6 = ra3*ra4+1.0;
|
---|
2730 | ra7 = 0.3149802624737183*pow(ra6,1.6666666666666667)+
|
---|
2731 | 0.3149802624737183*pow(ra5,1.6666666666666667);
|
---|
2732 | ra8 = 1/pow(ra7,0.5);
|
---|
2733 | ra9 = 1/(ra1*ra1);
|
---|
2734 | ra10 = 1/pow(ra1,1.6666666666666667);
|
---|
2735 | ra11 = 1/pow(ra1,0.66666666666666663);
|
---|
2736 | ra12 = 1/pow(ra1,0.33333333333333331);
|
---|
2737 | ra13 = 0.62035049089940009*c_3*ra12+0.38483473155912662*c_4*ra11+c_2;
|
---|
2738 | ra14 = 0.62035049089940009*c_5*ra12+0.38483473155912662*c_6*ra11
|
---|
2739 | +0.238732414637843*c_7*ra4+1.0;
|
---|
2740 | ra15 = 1/ra14;
|
---|
2741 | ra16 = (-(0.20678349696646667*c_3*ra2)-(0.25655648770608441*c_4*ra10))*ra15-
|
---|
2742 | ((1.0*(-(0.20678349696646667*c_5*ra2)-(0.25655648770608441*c_6*ra10)-
|
---|
2743 | (0.238732414637843*c_7*ra9))*ra13)/pow(ra14,2.0));
|
---|
2744 | ra17 = 1/pow(ra1,1.1666666666666667);
|
---|
2745 | ra18 = ra13*ra15+c_1;
|
---|
2746 | ra19 = 1/ra18;
|
---|
2747 | ra20 = exp(-1.0*a_*c_infin*mdr*ra17*ra19);
|
---|
2748 | double dp86c_drho_a = 0.79370052598409979*ra0*ra2*ra8*ra18*
|
---|
2749 | ((1.1666666666666667*a_*c_infin*mdr*ra19)/pow(ra1,2.1666666666666665)
|
---|
2750 | +(a_*c_infin*mdr*ra17*ra16)/ra18*ra18)*
|
---|
2751 | ra20-((1.0582673679787997*ra0*ra8*ra18*ra20)/
|
---|
2752 | pow(ra1,2.3333333333333335))-
|
---|
2753 | ((0.3968502629920499*ra0*ra2*(0.52496710412286385*(ra4-(1.0*ra3*ra9))*
|
---|
2754 | pow(ra6,0.66666666666666663)+0.52496710412286385*(-(1.0*ra4)+ra3*ra9)
|
---|
2755 | *pow(ra5,0.66666666666666663))*ra18*ra20)/pow(ra7,1.5))+
|
---|
2756 | 0.79370052598409979*ra0*ra2*ra8*ra16*ra20;
|
---|
2757 |
|
---|
2758 | return dp86c_drho_a;
|
---|
2759 | }
|
---|
2760 |
|
---|
2761 | double
|
---|
2762 | NewP86CFunctional::gab_deriv(double rho_a, double rho_b, double mdr)
|
---|
2763 | {
|
---|
2764 | double c_infin, c_1, c_2, c_3, c_4, c_5, c_6, c_7;
|
---|
2765 |
|
---|
2766 | c_infin = C1_ + C2_;
|
---|
2767 | c_1 = C1_; c_2 = C2_; c_3 = C3_; c_4 = C4_; c_5 = C5_; c_6 = C6_; c_7 = C7_;
|
---|
2768 |
|
---|
2769 | double g0,g1,g2,g3,g4,g5,g6,g7;
|
---|
2770 | g0 = rho_b+rho_a;
|
---|
2771 | g1 = -(1.0*rho_b)+rho_a;
|
---|
2772 | g2 = 1/g0;
|
---|
2773 | g3 = 1/pow(0.3149802624737183*pow(g1*g2+1.0,1.6666666666666667
|
---|
2774 | )+0.3149802624737183*pow(-(1.0*g1*g2)+1.0,
|
---|
2775 | 1.6666666666666667),0.5);
|
---|
2776 | g4 = 1/pow(g0,0.66666666666666663);
|
---|
2777 | g5 = 1/pow(g0,0.33333333333333331);
|
---|
2778 | g6 = (0.62035049089940009*c_3*g5+0.38483473155912662*c_4*g4+c_2)
|
---|
2779 | /pow(0.62035049089940009*c_5*g5+0.38483473155912662*c_6*g4+
|
---|
2780 | 0.238732414637843*c_7*g2+1.0,1.0)+c_1;
|
---|
2781 | g7 = exp(-1.0*a_*c_infin*mdr*1.0/pow(g0,1.1666666666666667)*1.0/g6);
|
---|
2782 | double dp86c_dmdr = (1.5874010519681996*mdr*g3*g6*g7)/pow(g0,1.3333333333333333)
|
---|
2783 | -((0.79370052598409979*a_*c_infin*pow(mdr,2.0)*g3*g7)/pow(g0,2.5));
|
---|
2784 |
|
---|
2785 | return dp86c_dmdr/mdr;
|
---|
2786 | }
|
---|
2787 |
|
---|
2788 | void
|
---|
2789 | NewP86CFunctional::point(const PointInputData &id,
|
---|
2790 | PointOutputData &od)
|
---|
2791 | {
|
---|
2792 | od.zero();
|
---|
2793 |
|
---|
2794 | // Precalculate terms for efficiency
|
---|
2795 | double rho = id.a.rho + id.b.rho;
|
---|
2796 | double rs = pow( (3./(4.*M_PI*rho)), (1./3.));
|
---|
2797 | double zet = (id.a.rho - id.b.rho)/rho;
|
---|
2798 | double c_infin = C1_ + C2_;
|
---|
2799 | double gamma_aa = id.a.gamma;
|
---|
2800 | double gamma_bb = id.b.gamma;
|
---|
2801 | double gamma_ab = id.gamma_ab;
|
---|
2802 | double mdr = sqrt(gamma_aa + gamma_bb + 2.*gamma_ab);
|
---|
2803 |
|
---|
2804 | double e0 = 1/pow(rho,0.66666666666666663);
|
---|
2805 | double e1 = 1/pow(rho,0.33333333333333331);
|
---|
2806 | double e2 = (0.62035049089940009*C3_*e1+0.38483473155912662*C4_*e0+C2_)
|
---|
2807 | /(0.62035049089940009*C5_*e1+0.38483473155912662*C6_*e0+(
|
---|
2808 | 0.238732414637843*C7_)/rho+1.0)+C1_;
|
---|
2809 | double fp86 = (0.79370052598409979*pow(mdr,2.0)*e2)*exp(-1.0*a_*
|
---|
2810 | c_infin*mdr*1.0/e2*1.0/pow(rho,1.1666666666666667)
|
---|
2811 | )/pow(rho,1.3333333333333333)/pow(0.3149802624737183*pow(zet+
|
---|
2812 | 1.0,1.6666666666666667)+0.3149802624737183*pow(-(1.0*zet)+
|
---|
2813 | 1.0,1.6666666666666667),0.5);
|
---|
2814 | od.energy += fp86;
|
---|
2815 |
|
---|
2816 | if (compute_potential_) {
|
---|
2817 | double dfp86_drhoa = rho_deriv(id.a.rho, id.b.rho, mdr);
|
---|
2818 | double dfp86_drhob = rho_deriv(id.b.rho, id.a.rho, mdr);
|
---|
2819 | double dfp86_dgab = gab_deriv(id.a.rho, id.b.rho, mdr);
|
---|
2820 |
|
---|
2821 | od.df_drho_a += dfp86_drhoa;
|
---|
2822 | od.df_drho_b += dfp86_drhob;
|
---|
2823 | od.df_dgamma_ab += dfp86_dgab;
|
---|
2824 | od.df_dgamma_aa += 0.5*od.df_dgamma_ab;
|
---|
2825 | od.df_dgamma_bb += 0.5*od.df_dgamma_ab;
|
---|
2826 | }
|
---|
2827 |
|
---|
2828 | }
|
---|
2829 |
|
---|
2830 | /////////////////////////////////////////////////////////////////////////////
|
---|
2831 | // PBECFunctional
|
---|
2832 |
|
---|
2833 | static ClassDesc PBECFunctional_cd(
|
---|
2834 | typeid(PBECFunctional),"PBECFunctional",1,"public DenFunctional",
|
---|
2835 | 0, create<PBECFunctional>, create<PBECFunctional>);
|
---|
2836 |
|
---|
2837 | PBECFunctional::PBECFunctional(StateIn& s):
|
---|
2838 | SavableState(s),
|
---|
2839 | DenFunctional(s)
|
---|
2840 | {
|
---|
2841 | local_ << SavableState::restore_state(s);
|
---|
2842 | s.get(gamma);
|
---|
2843 | s.get(beta);
|
---|
2844 | }
|
---|
2845 |
|
---|
2846 | PBECFunctional::PBECFunctional()
|
---|
2847 | {
|
---|
2848 | local_ = new PW92LCFunctional;
|
---|
2849 | local_->set_compute_potential(1);
|
---|
2850 |
|
---|
2851 | init_constants();
|
---|
2852 | }
|
---|
2853 |
|
---|
2854 | PBECFunctional::PBECFunctional(const Ref<KeyVal>& keyval):
|
---|
2855 | DenFunctional(keyval)
|
---|
2856 | {
|
---|
2857 | local_ << keyval->describedclassvalue("local");
|
---|
2858 | if (local_.null()) local_ = new PW92LCFunctional;
|
---|
2859 | local_->set_compute_potential(1);
|
---|
2860 |
|
---|
2861 | init_constants();
|
---|
2862 | gamma = keyval->doublevalue("gamma", KeyValValuedouble(gamma));
|
---|
2863 | beta = keyval->doublevalue("beta", KeyValValuedouble(beta));
|
---|
2864 | }
|
---|
2865 |
|
---|
2866 | void
|
---|
2867 | PBECFunctional::init_constants()
|
---|
2868 | {
|
---|
2869 | // in paper
|
---|
2870 | // gamma = 0.031091
|
---|
2871 | // beta = 0.066725
|
---|
2872 | // in PBE.f:
|
---|
2873 | gamma = 0.03109069086965489503494086371273;
|
---|
2874 | beta = 0.06672455060314922;
|
---|
2875 | }
|
---|
2876 |
|
---|
2877 | PBECFunctional::~PBECFunctional()
|
---|
2878 | {
|
---|
2879 | }
|
---|
2880 |
|
---|
2881 | void
|
---|
2882 | PBECFunctional::save_data_state(StateOut& s)
|
---|
2883 | {
|
---|
2884 | DenFunctional::save_data_state(s);
|
---|
2885 | SavableState::save_state(local_.pointer(),s);
|
---|
2886 | s.put(gamma);
|
---|
2887 | s.put(beta);
|
---|
2888 | }
|
---|
2889 |
|
---|
2890 | int
|
---|
2891 | PBECFunctional::need_density_gradient()
|
---|
2892 | {
|
---|
2893 | return 1;
|
---|
2894 | }
|
---|
2895 |
|
---|
2896 | void
|
---|
2897 | PBECFunctional::set_spin_polarized(int a)
|
---|
2898 | {
|
---|
2899 | spin_polarized_ = a;
|
---|
2900 | local_->set_spin_polarized(a);
|
---|
2901 | }
|
---|
2902 |
|
---|
2903 | double
|
---|
2904 | PBECFunctional::rho_deriv(double rho_a, double rho_b, double mdr,
|
---|
2905 | double ec_local, double ec_local_dra)
|
---|
2906 | {
|
---|
2907 | if (mdr < MIN_SQRTGAMMA) {
|
---|
2908 | return 0;
|
---|
2909 | }
|
---|
2910 |
|
---|
2911 | if (rho_b < MIN_DENSITY) {
|
---|
2912 | #define Log(x) log(x)
|
---|
2913 | #define Power(x,y) pow(x,y)
|
---|
2914 | #define Pi M_PI
|
---|
2915 | #define E M_E
|
---|
2916 | double ec = ec_local;
|
---|
2917 | double decdrhoa = ec_local_dra;
|
---|
2918 | double rhoa = rho_a;
|
---|
2919 | double result =
|
---|
2920 | (gamma*((-2*beta*Power(mdr,2)*Power(Pi,0.3333333333333333)*rhoa*
|
---|
2921 | (Power(beta,2)*decdrhoa*Power(E,(4*ec)/gamma)*Power(mdr,4)*
|
---|
2922 | Power(Pi,0.6666666666666666)*
|
---|
2923 | (Power(6,0.6666666666666666)*beta*Power(E,(2*ec)/gamma)*
|
---|
2924 | Power(mdr,2)*Power(Pi,0.3333333333333333) -
|
---|
2925 | 96*(-1 + Power(E,(2*ec)/gamma))*gamma*
|
---|
2926 | Power(rhoa,2.3333333333333335)) +
|
---|
2927 | 896*Power(6,0.3333333333333333)*
|
---|
2928 | Power(-1 + Power(E,(2*ec)/gamma),3)*Power(gamma,3)*
|
---|
2929 | Power(rhoa,3.6666666666666665)*
|
---|
2930 | (-(beta*Power(E,(2*ec)/gamma)*Power(mdr,2)*
|
---|
2931 | Power(Pi,0.3333333333333333)) +
|
---|
2932 | 4*Power(6,0.3333333333333333)*(-1 + Power(E,(2*ec)/gamma))*
|
---|
2933 | gamma*Power(rhoa,2.3333333333333335))))/
|
---|
2934 | (gamma*(Power(6,0.3333333333333333)*Power(beta,2)*
|
---|
2935 | Power(E,(2*ec)/gamma)*Power(mdr,4)*Power(Pi,0.6666666666666666)\
|
---|
2936 | - 8*Power(6,0.6666666666666666)*beta*
|
---|
2937 | (-1 + Power(E,(2*ec)/gamma))*gamma*Power(mdr,2)*
|
---|
2938 | Power(Pi,0.3333333333333333)*Power(rhoa,2.3333333333333335) +
|
---|
2939 | 384*Power(-1 + Power(E,(2*ec)/gamma),2)*Power(gamma,2)*
|
---|
2940 | Power(rhoa,4.666666666666667))*
|
---|
2941 | (Power(6,0.3333333333333333)*Power(beta,2)*Power(E,(4*ec)/gamma)*
|
---|
2942 | Power(mdr,4)*Power(Pi,0.6666666666666666) -
|
---|
2943 | 8*Power(6,0.6666666666666666)*beta*Power(E,(2*ec)/gamma)*
|
---|
2944 | (-1 + Power(E,(2*ec)/gamma))*gamma*Power(mdr,2)*
|
---|
2945 | Power(Pi,0.3333333333333333)*Power(rhoa,2.3333333333333335) +
|
---|
2946 | 384*Power(-1 + Power(E,(2*ec)/gamma),2)*Power(gamma,2)*
|
---|
2947 | Power(rhoa,4.666666666666667))) +
|
---|
2948 | Log(1 + (beta*(-1 + Power(E,(2*ec)/gamma))*Power(mdr,2)*
|
---|
2949 | Power(Pi/6.,0.3333333333333333)*
|
---|
2950 | (Power(6,0.6666666666666666)*beta*Power(E,(2*ec)/gamma)*
|
---|
2951 | Power(mdr,2)*Power(Pi,0.3333333333333333) -
|
---|
2952 | 48*(-1 + Power(E,(2*ec)/gamma))*gamma*
|
---|
2953 | Power(rhoa,2.3333333333333335)))/
|
---|
2954 | (-(Power(6,0.3333333333333333)*Power(beta,2)*Power(E,(4*ec)/gamma)*
|
---|
2955 | Power(mdr,4)*Power(Pi,0.6666666666666666)) +
|
---|
2956 | 8*Power(6,0.6666666666666666)*beta*Power(E,(2*ec)/gamma)*
|
---|
2957 | (-1 + Power(E,(2*ec)/gamma))*gamma*Power(mdr,2)*
|
---|
2958 | Power(Pi,0.3333333333333333)*Power(rhoa,2.3333333333333335) -
|
---|
2959 | 384*Power(-1 + Power(E,(2*ec)/gamma),2)*Power(gamma,2)*
|
---|
2960 | Power(rhoa,4.666666666666667)))))/2.;
|
---|
2961 | return result;
|
---|
2962 | }
|
---|
2963 | if (rho_a < MIN_DENSITY) {
|
---|
2964 | // df_drho_a diverges for this case
|
---|
2965 | return 0.0;
|
---|
2966 | }
|
---|
2967 |
|
---|
2968 | double ra0,ra1,ra2,ra3,ra4,ra5,ra6,ra7,ra8,ra9,ra10,ra11,ra12,ra13,ra14,
|
---|
2969 | ra15,ra16,ra17,ra18,ra19,ra20,ra21,ra22,ra23,ra24,ra25,ra26,ra27,ra28,
|
---|
2970 | ra29,ra30,ra31,ra32,ra33,ra34,ra35;
|
---|
2971 | double dpbec_drho_a;
|
---|
2972 |
|
---|
2973 | ra0 = rho_b+rho_a;
|
---|
2974 | ra1 = -(1.0*rho_b)+rho_a;
|
---|
2975 | ra2 = 1/pow(ra0,2.0);
|
---|
2976 | ra3 = 1/pow(ra0,1.0);
|
---|
2977 | ra4 = -(1.0*ra1*ra3)+1.0;
|
---|
2978 | ra5 = ra1*ra3+1.0;
|
---|
2979 | ra6 = (0.66666666666666663*(ra3-(1.0*ra1*ra2)))/pow(ra5,
|
---|
2980 | 0.33333333333333331)+(0.66666666666666663*(-(1.0*ra3)+ra1*ra2
|
---|
2981 | ))/pow(ra4,0.33333333333333331);
|
---|
2982 | ra7 = pow(ra5,0.66666666666666663)+pow(ra4,0.66666666666666663);
|
---|
2983 | ra8 = pow(ra7,2.0);
|
---|
2984 | ra9 = pow(mdr,2.0);
|
---|
2985 | ra10 = 1/pow(ra0,2.3333333333333335);
|
---|
2986 | ra11 = 1/ra8;
|
---|
2987 | ra12 = pow(ra7,3.0);
|
---|
2988 | ra13 = 1/ra12;
|
---|
2989 | ra14 = ec_local;
|
---|
2990 | ra15 = 1/pow(gamma,1.0);
|
---|
2991 | ra16 = 1/exp(8.0*ra13*ra14*ra15);
|
---|
2992 | ra17 = ra16-1.0;
|
---|
2993 | ra18 = 1/pow(ra17,1.0);
|
---|
2994 | ra19 = 0.25387282439081477*beta*ra9*ra10*ra11*ra18*ra15;
|
---|
2995 | ra20 = ra19+1.0;
|
---|
2996 | ra21 = pow(beta,2.0);
|
---|
2997 | ra22 = pow(mdr,4.0);
|
---|
2998 | ra23 = 1/pow(ra0,4.666666666666667);
|
---|
2999 | ra24 = 1/pow(ra7,4.0);
|
---|
3000 | ra25 = 1/pow(ra17,2.0);
|
---|
3001 | ra26 = 1/pow(gamma,2.0);
|
---|
3002 | ra27 = ra19+0.064451410964169495*ra21*ra22*ra23*ra24*ra25*ra26+
|
---|
3003 | 1.0;
|
---|
3004 | ra28 = 1/pow(ra27,1.0);
|
---|
3005 | ra29 = 0.25387282439081477*beta*ra9*ra10*ra11*ra20*ra28*ra15+1.0
|
---|
3006 | ;
|
---|
3007 | ra30 = log(ra29);
|
---|
3008 | ra31 = 1/pow(ra0,3.3333333333333335);
|
---|
3009 | ra32 = -(0.50774564878162953*beta*ra9*ra10*ra6*ra13*ra18*ra15);
|
---|
3010 | ra33 = -(0.59236992357856788*beta*ra9*ra31*ra11*ra18*ra15);
|
---|
3011 | ra34 = -(8.0*ra13*ec_local_dra*ra15)+24.0*ra6*
|
---|
3012 | ra24*ra14*ra15;
|
---|
3013 | ra35 = -(0.25387282439081477*beta*ra9*ra10*ra11*ra25*ra16*ra34*
|
---|
3014 | ra15);
|
---|
3015 | dpbec_drho_a = (0.125*ra0*ra12*(-((0.25387282439081477*beta*ra9*
|
---|
3016 | ra10*ra11*ra20*(ra35+ra33+ra32-((0.12890282192833899*ra21*ra22*
|
---|
3017 | ra23*ra24*ra16*ra34*ra26)/pow(ra17,3.0))-((0.30077325116612436*
|
---|
3018 | ra21*ra22*ra24*ra25*ra26)/pow(ra0,5.666666666666667))-((
|
---|
3019 | 0.25780564385667798*ra21*ra22*ra23*ra6*ra25*ra26)/pow(ra7,5.0))
|
---|
3020 | )*ra15)/pow(ra27,2.0))+0.25387282439081477*beta*ra9*ra10*ra11*
|
---|
3021 | ra28*(ra35+ra33+ra32)*ra15-(0.59236992357856788*beta*ra9*ra31*
|
---|
3022 | ra11*ra20*ra28*ra15)-(0.50774564878162953*beta*ra9*ra10*ra6*ra13*
|
---|
3023 | ra20*ra28*ra15))*gamma)/pow(ra29,1.0)+0.125*ra12*ra30*gamma+
|
---|
3024 | 0.375*ra0*ra6*ra8*ra30*gamma;
|
---|
3025 |
|
---|
3026 | return dpbec_drho_a;
|
---|
3027 | }
|
---|
3028 |
|
---|
3029 | double
|
---|
3030 | PBECFunctional::gab_deriv(double rho, double phi, double mdr, double ec_local)
|
---|
3031 | {
|
---|
3032 | if (rho < MIN_DENSITY) return 0;
|
---|
3033 |
|
---|
3034 | if (mdr < MIN_SQRTGAMMA) {
|
---|
3035 | double result
|
---|
3036 | = (beta*phi*pow(M_PI/3.,1./3.))/(8.*pow(rho,4./3.));
|
---|
3037 | return result;
|
---|
3038 | }
|
---|
3039 |
|
---|
3040 | double g0,g1,g2,g3,g4,g5,g6,g7,g8,g9,g10,g11,g12,g13,g14,g15,g16;
|
---|
3041 | double dpbec_dmdr;
|
---|
3042 | g0 = pow(phi,3.0);
|
---|
3043 | g1 = pow(beta,2.0);
|
---|
3044 | g2 = pow(mdr,3.0);
|
---|
3045 | g3 = 1/pow(phi,4.0);
|
---|
3046 | g4 = 1/pow(rho,4.666666666666667);
|
---|
3047 | g5 = 1/pow(gamma,1.0);
|
---|
3048 | g6 = 1/exp(1.0*1.0/g0*ec_local*g5)-1.0;
|
---|
3049 | g7 = 1/pow(g6,1.0);
|
---|
3050 | g8 = 1/pow(g6,2.0);
|
---|
3051 | g9 = 1/pow(gamma,2.0);
|
---|
3052 | g10 = pow(mdr,2.0);
|
---|
3053 | g11 = 1/pow(phi,2.0);
|
---|
3054 | g12 = 1/pow(rho,2.3333333333333335);
|
---|
3055 | g13 = 0.063468206097703692*beta*g10*g11*g12*g7*g5;
|
---|
3056 | g14 = g13+0.0040282131852605934*g1*pow(mdr,4.0)*g3*g4*g8*g9+
|
---|
3057 | 1.0;
|
---|
3058 | g15 = 1/pow(g14,1.0);
|
---|
3059 | g16 = g13+1.0;
|
---|
3060 | dpbec_dmdr = (g0*rho*(0.12693641219540738*beta*mdr*g11*g12*g16*g15
|
---|
3061 | *g5-((0.063468206097703692*beta*g10*g11*g12*(
|
---|
3062 | 0.12693641219540738*beta*mdr*g11*g12*g7*g5+0.016112852741042374
|
---|
3063 | *g1*g2*g3*g4*g8*g9)*g16*g5)/pow(g14,2.0))+0.0080564263705211869
|
---|
3064 | *g1*g2*g3*g4*g7*g15*g9)*gamma)/pow(0.063468206097703692*beta*g10*
|
---|
3065 | g11*g12*g16*g15*g5+1.0,1.0);
|
---|
3066 |
|
---|
3067 | return dpbec_dmdr/mdr;
|
---|
3068 | }
|
---|
3069 |
|
---|
3070 | void
|
---|
3071 | PBECFunctional::point(const PointInputData &id,
|
---|
3072 | PointOutputData &od)
|
---|
3073 | {
|
---|
3074 | double ec_local, dec_local_rs, dec_local_zeta;
|
---|
3075 | local_->point_lc(id, od, ec_local, dec_local_rs, dec_local_zeta);
|
---|
3076 | double rho = id.a.rho+id.b.rho;
|
---|
3077 | double rho_13 = pow(rho, 1./3.);
|
---|
3078 | double rho_43 = rho*rho_13;
|
---|
3079 | double zeta = (id.a.rho - id.b.rho)/rho;
|
---|
3080 | double phi = 0.5*(pow(1+zeta, 2./3.)+pow(1-zeta, 2./3.));
|
---|
3081 | double mdr = sqrt(id.a.gamma + id.b.gamma + 2*id.gamma_ab);
|
---|
3082 |
|
---|
3083 | double pbec;
|
---|
3084 |
|
---|
3085 | double e0,e1,e2,e3,e4,e5,e6;
|
---|
3086 | e0 = pow(phi,3.0);
|
---|
3087 | e1 = pow(mdr,2.0);
|
---|
3088 | e2 = 1/pow(phi,2.0);
|
---|
3089 | e3 = 1/pow(rho,2.3333333333333335);
|
---|
3090 | e4 = 1/pow(gamma,1.0);
|
---|
3091 | e5 = 1/exp(1.0*1.0/e0*ec_local*e4)-1.0;
|
---|
3092 | e6 = (0.063468206097703692*beta*e1*e2*e3*e4)/pow(e5,1.0);
|
---|
3093 | pbec = e0*rho*log((0.063468206097703692*beta*e1*e2*e3*(e6+1.0)*
|
---|
3094 | e4)/pow(e6+(0.0040282131852605934*pow(beta,2.0)*pow(mdr,4.0))
|
---|
3095 | /pow(phi,4.0)/pow(rho,4.666666666666667)/pow(e5,2.0)/pow(
|
---|
3096 | gamma,2.0)+1.0,1.0)+1.0)*gamma;
|
---|
3097 |
|
---|
3098 | if (compute_potential_) {
|
---|
3099 | double drs_drho_a = -0.20678349696646667/rho_43; // == drs_drho_b
|
---|
3100 | double dzeta_drho_a = 0.;
|
---|
3101 | if (zeta < MAX_ZETA) dzeta_drho_a = 1./rho * ( 1. - zeta);
|
---|
3102 | double dzeta_drho_b = 0.;
|
---|
3103 | if (zeta > MIN_ZETA) dzeta_drho_b = 1./rho * (-1. - zeta);
|
---|
3104 |
|
---|
3105 | //double ec_local_dra = od.df_drho_a;
|
---|
3106 | //double ec_local_drb = od.df_drho_b;
|
---|
3107 | double ec_local_dra
|
---|
3108 | = dec_local_rs*drs_drho_a + dec_local_zeta*dzeta_drho_a;
|
---|
3109 | double ec_local_drb
|
---|
3110 | = dec_local_rs*drs_drho_a + dec_local_zeta*dzeta_drho_b;
|
---|
3111 |
|
---|
3112 | double df_drho_a = rho_deriv(id.a.rho, id.b.rho, mdr,
|
---|
3113 | ec_local, ec_local_dra);
|
---|
3114 | double df_drho_b = rho_deriv(id.b.rho, id.a.rho, mdr,
|
---|
3115 | ec_local, ec_local_drb);
|
---|
3116 |
|
---|
3117 | double df_dgab = gab_deriv(rho, phi, mdr, ec_local);
|
---|
3118 |
|
---|
3119 | od.df_drho_a += df_drho_a;
|
---|
3120 | od.df_drho_b += df_drho_b;
|
---|
3121 | od.df_dgamma_aa += 0.5*df_dgab;
|
---|
3122 | od.df_dgamma_bb += 0.5*df_dgab;
|
---|
3123 | od.df_dgamma_ab += df_dgab;
|
---|
3124 | }
|
---|
3125 |
|
---|
3126 | od.energy += pbec;
|
---|
3127 |
|
---|
3128 | // cout << scprintf("id.a.del_rho = %12.8f %12.8f %12.8f", id.a.del_rho[0],
|
---|
3129 | // id.a.del_rho[1], id.a.del_rho[2])
|
---|
3130 | // << endl;
|
---|
3131 | // cout << scprintf("id.b.del_rho = %12.8f %12.8f %12.8f", id.b.del_rho[0],
|
---|
3132 | // id.b.del_rho[1], id.b.del_rho[2])
|
---|
3133 | // << endl;
|
---|
3134 | // cout << "id.a.gamma = " << id.a.gamma << endl;
|
---|
3135 | // cout << "id.b.gamma = " << id.b.gamma << endl;
|
---|
3136 | // cout << "od.df_drho_a = " << od.df_drho_a << endl;
|
---|
3137 | // cout << "od.df_drho_b = " << od.df_drho_b << endl;
|
---|
3138 | // cout << "od.df_dgamma_aa = " << od.df_dgamma_aa << endl;
|
---|
3139 | // cout << "od.df_dgamma_ab = " << od.df_dgamma_ab << endl;
|
---|
3140 | // cout << "od.df_dgamma_bb = " << od.df_dgamma_bb << endl;
|
---|
3141 | }
|
---|
3142 |
|
---|
3143 | /////////////////////////////////////////////////////////////////////////////
|
---|
3144 | // PW91CFunctional
|
---|
3145 |
|
---|
3146 | static ClassDesc PW91CFunctional_cd(
|
---|
3147 | typeid(PW91CFunctional),"PW91CFunctional",1,"public DenFunctional",
|
---|
3148 | 0, create<PW91CFunctional>, create<PW91CFunctional>);
|
---|
3149 |
|
---|
3150 | PW91CFunctional::PW91CFunctional(StateIn& s):
|
---|
3151 | SavableState(s),
|
---|
3152 | DenFunctional(s)
|
---|
3153 | {
|
---|
3154 | local_ << SavableState::restore_state(s);
|
---|
3155 | init_constants();
|
---|
3156 | }
|
---|
3157 |
|
---|
3158 | PW91CFunctional::PW91CFunctional()
|
---|
3159 | {
|
---|
3160 | local_ = new PW92LCFunctional;
|
---|
3161 | local_->set_compute_potential(1);
|
---|
3162 | init_constants();
|
---|
3163 | }
|
---|
3164 |
|
---|
3165 | PW91CFunctional::PW91CFunctional(const Ref<KeyVal>& keyval):
|
---|
3166 | DenFunctional(keyval)
|
---|
3167 | {
|
---|
3168 | local_ << keyval->describedclassvalue("local");
|
---|
3169 | if (local_.null()) local_ = new PW92LCFunctional;
|
---|
3170 | local_->set_compute_potential(1);
|
---|
3171 | init_constants();
|
---|
3172 | }
|
---|
3173 |
|
---|
3174 | void
|
---|
3175 | PW91CFunctional::init_constants()
|
---|
3176 | {
|
---|
3177 | a = 23.266;
|
---|
3178 | b = 7.389e-3;
|
---|
3179 | c = 8.723;
|
---|
3180 | d = 0.472;
|
---|
3181 | alpha = 0.09;
|
---|
3182 | c_c0 = 0.004235;
|
---|
3183 | // c_x = -0.001667 in Phys Rev B 46, p 6671, Perdew, et. al.
|
---|
3184 | // c_x as in PBE.f:
|
---|
3185 | c_x = -0.001667212;
|
---|
3186 | }
|
---|
3187 |
|
---|
3188 | PW91CFunctional::~PW91CFunctional()
|
---|
3189 | {
|
---|
3190 | }
|
---|
3191 |
|
---|
3192 | void
|
---|
3193 | PW91CFunctional::save_data_state(StateOut& s)
|
---|
3194 | {
|
---|
3195 | DenFunctional::save_data_state(s);
|
---|
3196 | SavableState::save_state(local_.pointer(),s);
|
---|
3197 | }
|
---|
3198 |
|
---|
3199 | int
|
---|
3200 | PW91CFunctional::need_density_gradient()
|
---|
3201 | {
|
---|
3202 | return 1;
|
---|
3203 | }
|
---|
3204 |
|
---|
3205 | void
|
---|
3206 | PW91CFunctional::set_spin_polarized(int a)
|
---|
3207 | {
|
---|
3208 | spin_polarized_ = a;
|
---|
3209 | local_->set_spin_polarized(a);
|
---|
3210 | }
|
---|
3211 |
|
---|
3212 | double
|
---|
3213 | PW91CFunctional::limit_df_drhoa(double rhoa, double mdr,
|
---|
3214 | double ec, double decdrhoa)
|
---|
3215 | {
|
---|
3216 | double result;
|
---|
3217 | // blame mathematica
|
---|
3218 | double v = 15.755920349483144;
|
---|
3219 | double cx = c_x;
|
---|
3220 | double cc0 = c_c0;
|
---|
3221 | double beta = v * cc0;
|
---|
3222 | double e2ec = Power(E,(2*alpha*ec)/Power(beta,2));
|
---|
3223 | double e4ec = e2ec * e2ec;
|
---|
3224 | double e8ec = e4ec * e4ec;
|
---|
3225 | double e25mdr2 = Power(E,-(25*Power(mdr,2))/
|
---|
3226 | (Power(6,0.6666666666666666)*Power(Pi,1.3333333333333333)*
|
---|
3227 | Power(rhoa,2.6666666666666665)));
|
---|
3228 | result =
|
---|
3229 | (Power(mdr,2)*Power(Pi/6.,0.3333333333333333)*
|
---|
3230 | (1750*Power(6,0.3333333333333333)*a*Power(Pi,0.6666666666666666) -
|
---|
3231 | 1750000*Power(6,0.3333333333333333)*c*cc0*
|
---|
3232 | Power(Pi,0.6666666666666666) -
|
---|
3233 | 2500000*Power(6,0.3333333333333333)*c*cx*
|
---|
3234 | Power(Pi,0.6666666666666666) +
|
---|
3235 | (8988*Pi)/Power(1/rhoa,0.3333333333333333) -
|
---|
3236 | (3500000*cc0*Pi)/Power(1/rhoa,0.3333333333333333) -
|
---|
3237 | (5000000*cx*Pi)/Power(1/rhoa,0.3333333333333333) +
|
---|
3238 | 875*Power(6,0.6666666666666666)*b*Power(Pi,0.3333333333333333)*
|
---|
3239 | Power(1/rhoa,0.3333333333333333) -
|
---|
3240 | 875000*Power(6,0.6666666666666666)*cc0*d*
|
---|
3241 | Power(Pi,0.3333333333333333)*Power(1/rhoa,0.3333333333333333) -
|
---|
3242 | 1250000*Power(6,0.6666666666666666)*cx*d*
|
---|
3243 | Power(Pi,0.3333333333333333)*Power(1/rhoa,0.3333333333333333) -
|
---|
3244 | 26250000*b*cc0*Power(1/rhoa,0.6666666666666666) -
|
---|
3245 | 37500000*b*cx*Power(1/rhoa,0.6666666666666666))*v)
|
---|
3246 | *e25mdr2/(1.4e7*
|
---|
3247 | (2*Power(6,0.3333333333333333)*c*Power(Pi,0.6666666666666666) +
|
---|
3248 | (4*Pi)/Power(1/rhoa,0.3333333333333333) +
|
---|
3249 | Power(6,0.6666666666666666)*d*Power(Pi,0.3333333333333333)*
|
---|
3250 | Power(1/rhoa,0.3333333333333333) +
|
---|
3251 | 30*b*Power(1/rhoa,0.6666666666666666))*Power(rhoa,2.3333333333333335)
|
---|
3252 | );
|
---|
3253 | result +=
|
---|
3254 | rhoa*((Power(beta,2)*
|
---|
3255 | ((-2*alpha*(-(alpha*e4ec*Power(mdr,4)*
|
---|
3256 | Power(Pi/6.,0.6666666666666666))/
|
---|
3257 | (32.*beta*(-1 + e4ec)*Power(rhoa,4.666666666666667)) +
|
---|
3258 | (Power(mdr,2)*Power(Pi/6.,0.3333333333333333))/
|
---|
3259 | (8.*Power(rhoa,2.3333333333333335)))*
|
---|
3260 | ((-7*Power(alpha,2)*e8ec*Power(mdr,4)*
|
---|
3261 | Power(Pi/6.,0.6666666666666666))/
|
---|
3262 | (24.*beta*Power(-1 + e4ec,2)*Power(rhoa,5.666666666666667))\
|
---|
3263 | - (Power(alpha,3)*decdrhoa*e8ec*Power(mdr,4)*
|
---|
3264 | Power(Pi/6.,0.6666666666666666))/
|
---|
3265 | (2.*Power(beta,3)*Power(-1 + e4ec,3)*
|
---|
3266 | Power(rhoa,4.666666666666667)) +
|
---|
3267 | (7*alpha*e4ec*Power(mdr,2)*Power(Pi/6.,0.3333333333333333))/
|
---|
3268 | (12.*(-1 + e4ec)*Power(rhoa,3.3333333333333335)) +
|
---|
3269 | (Power(alpha,2)*decdrhoa*e4ec*Power(mdr,2)*
|
---|
3270 | Power(Pi/6.,0.3333333333333333))/
|
---|
3271 | (Power(beta,2)*Power(-1 + e4ec,2)*
|
---|
3272 | Power(rhoa,2.3333333333333335))))/
|
---|
3273 | Power(beta + (Power(alpha,2)*e8ec*Power(mdr,4)*
|
---|
3274 | Power(Pi/6.,0.6666666666666666))/
|
---|
3275 | (16.*beta*Power(-1 + e4ec,2)*Power(rhoa,4.666666666666667)) -
|
---|
3276 | (alpha*e4ec*Power(mdr,2)*Power(Pi/6.,0.3333333333333333))/
|
---|
3277 | (4.*(-1 + e4ec)*Power(rhoa,2.3333333333333335)),2) +
|
---|
3278 | (2*alpha*((7*alpha*e4ec*Power(mdr,4)*
|
---|
3279 | Power(Pi/6.,0.6666666666666666))/
|
---|
3280 | (48.*beta*(-1 + e4ec)*Power(rhoa,5.666666666666667)) +
|
---|
3281 | (Power(alpha,2)*decdrhoa*e4ec*Power(mdr,4)*
|
---|
3282 | Power(Pi/6.,0.6666666666666666))/
|
---|
3283 | (8.*Power(beta,3)*Power(-1 + e4ec,2)*
|
---|
3284 | Power(rhoa,4.666666666666667)) -
|
---|
3285 | (7*Power(mdr,2)*Power(Pi/6.,0.3333333333333333))/
|
---|
3286 | (24.*Power(rhoa,3.3333333333333335))))/
|
---|
3287 | (beta + (Power(alpha,2)*e8ec*Power(mdr,4)*
|
---|
3288 | Power(Pi/6.,0.6666666666666666))/
|
---|
3289 | (16.*beta*Power(-1 + e4ec,2)*Power(rhoa,4.666666666666667)) -
|
---|
3290 | (alpha*e4ec*Power(mdr,2)*Power(Pi/6.,0.3333333333333333))/
|
---|
3291 | (4.*(-1 + e4ec)*Power(rhoa,2.3333333333333335)))))/
|
---|
3292 | (4.*alpha*(1 + (2*alpha*
|
---|
3293 | (-(alpha*e4ec*Power(mdr,4)*Power(Pi/6.,0.6666666666666666))/
|
---|
3294 | (32.*beta*(-1 + e4ec)*Power(rhoa,4.666666666666667)) +
|
---|
3295 | (Power(mdr,2)*Power(Pi/6.,0.3333333333333333))/
|
---|
3296 | (8.*Power(rhoa,2.3333333333333335))))/
|
---|
3297 | (beta + (Power(alpha,2)*e8ec*Power(mdr,4)*
|
---|
3298 | Power(Pi/6.,0.6666666666666666))/
|
---|
3299 | (16.*beta*Power(-1 + e4ec,2)*Power(rhoa,4.666666666666667)) -
|
---|
3300 | (alpha*e4ec*Power(mdr,2)*Power(Pi/6.,0.3333333333333333))/
|
---|
3301 | (4.*(-1 + e4ec)*Power(rhoa,2.3333333333333335))))) +
|
---|
3302 | (Power(mdr,4)*(1750*Power(6,0.3333333333333333)*a*
|
---|
3303 | Power(Pi,0.6666666666666666) -
|
---|
3304 | 1750000*Power(6,0.3333333333333333)*c*cc0*
|
---|
3305 | Power(Pi,0.6666666666666666) -
|
---|
3306 | 2500000*Power(6,0.3333333333333333)*c*cx*
|
---|
3307 | Power(Pi,0.6666666666666666) +
|
---|
3308 | (8988*Pi)/Power(1/rhoa,0.3333333333333333) -
|
---|
3309 | (3500000*cc0*Pi)/Power(1/rhoa,0.3333333333333333) -
|
---|
3310 | (5000000*cx*Pi)/Power(1/rhoa,0.3333333333333333) +
|
---|
3311 | 875*Power(6,0.6666666666666666)*b*Power(Pi,0.3333333333333333)*
|
---|
3312 | Power(1/rhoa,0.3333333333333333) -
|
---|
3313 | 875000*Power(6,0.6666666666666666)*cc0*d*
|
---|
3314 | Power(Pi,0.3333333333333333)*Power(1/rhoa,0.3333333333333333) -
|
---|
3315 | 1250000*Power(6,0.6666666666666666)*cx*d*
|
---|
3316 | Power(Pi,0.3333333333333333)*Power(1/rhoa,0.3333333333333333) -
|
---|
3317 | 26250000*b*cc0*Power(1/rhoa,0.6666666666666666) -
|
---|
3318 | 37500000*b*cx*Power(1/rhoa,0.6666666666666666))*v)
|
---|
3319 | *e25mdr2/(1.26e6*Pi*
|
---|
3320 | (2*Power(6,0.3333333333333333)*c*Power(Pi,0.6666666666666666) +
|
---|
3321 | (4*Pi)/Power(1/rhoa,0.3333333333333333) +
|
---|
3322 | Power(6,0.6666666666666666)*d*Power(Pi,0.3333333333333333)*
|
---|
3323 | Power(1/rhoa,0.3333333333333333) +
|
---|
3324 | 30*b*Power(1/rhoa,0.6666666666666666))*Power(rhoa,6)) -
|
---|
3325 | (Power(mdr,2)*Power(Pi/6.,0.3333333333333333)*
|
---|
3326 | (1750*Power(6,0.3333333333333333)*a*Power(Pi,0.6666666666666666) -
|
---|
3327 | 1750000*Power(6,0.3333333333333333)*c*cc0*
|
---|
3328 | Power(Pi,0.6666666666666666) -
|
---|
3329 | 2500000*Power(6,0.3333333333333333)*c*cx*
|
---|
3330 | Power(Pi,0.6666666666666666) +
|
---|
3331 | (8988*Pi)/Power(1/rhoa,0.3333333333333333) -
|
---|
3332 | (3500000*cc0*Pi)/Power(1/rhoa,0.3333333333333333) -
|
---|
3333 | (5000000*cx*Pi)/Power(1/rhoa,0.3333333333333333) +
|
---|
3334 | 875*Power(6,0.6666666666666666)*b*Power(Pi,0.3333333333333333)*
|
---|
3335 | Power(1/rhoa,0.3333333333333333) -
|
---|
3336 | 875000*Power(6,0.6666666666666666)*cc0*d*
|
---|
3337 | Power(Pi,0.3333333333333333)*Power(1/rhoa,0.3333333333333333) -
|
---|
3338 | 1250000*Power(6,0.6666666666666666)*cx*d*
|
---|
3339 | Power(Pi,0.3333333333333333)*Power(1/rhoa,0.3333333333333333) -
|
---|
3340 | 26250000*b*cc0*Power(1/rhoa,0.6666666666666666) -
|
---|
3341 | 37500000*b*cx*Power(1/rhoa,0.6666666666666666))*v)*e25mdr2/
|
---|
3342 | (6.e6*
|
---|
3343 | (2*Power(6,0.3333333333333333)*c*Power(Pi,0.6666666666666666) +
|
---|
3344 | (4*Pi)/Power(1/rhoa,0.3333333333333333) +
|
---|
3345 | Power(6,0.6666666666666666)*d*Power(Pi,0.3333333333333333)*
|
---|
3346 | Power(1/rhoa,0.3333333333333333) +
|
---|
3347 | 30*b*Power(1/rhoa,0.6666666666666666))*
|
---|
3348 | Power(rhoa,3.3333333333333335)) +
|
---|
3349 | (Power(mdr,2)*Power(Pi/6.,0.3333333333333333)*
|
---|
3350 | (1875*Power(6,0.6666666666666666)*Power(b,2)*
|
---|
3351 | Power(Pi,0.3333333333333333) -
|
---|
3352 | (500*Power(6,0.3333333333333333)*a*Power(Pi,1.6666666666666667))/
|
---|
3353 | Power(1/rhoa,1.3333333333333333) +
|
---|
3354 | (1284*Power(6,0.3333333333333333)*c*Power(Pi,1.6666666666666667))/
|
---|
3355 | Power(1/rhoa,1.3333333333333333) +
|
---|
3356 | (57780*b*Pi)/Power(1/rhoa,0.6666666666666666) -
|
---|
3357 | (750*b*c*Pi)/Power(1/rhoa,0.6666666666666666) +
|
---|
3358 | (750*a*d*Pi)/Power(1/rhoa,0.6666666666666666) +
|
---|
3359 | (7500*Power(6,0.3333333333333333)*a*b*
|
---|
3360 | Power(Pi,0.6666666666666666))/Power(1/rhoa,0.3333333333333333)\
|
---|
3361 | - 500*Power(6,0.6666666666666666)*b*Power(Pi,1.3333333333333333)*
|
---|
3362 | rhoa + 1284*Power(6,0.6666666666666666)*d*
|
---|
3363 | Power(Pi,1.3333333333333333)*rhoa)*v)*e25mdr2/
|
---|
3364 | (3.e6*
|
---|
3365 | Power(2*Power(6,0.3333333333333333)*c*
|
---|
3366 | Power(Pi,0.6666666666666666) +
|
---|
3367 | (4*Pi)/Power(1/rhoa,0.3333333333333333) +
|
---|
3368 | Power(6,0.6666666666666666)*d*Power(Pi,0.3333333333333333)*
|
---|
3369 | Power(1/rhoa,0.3333333333333333) +
|
---|
3370 | 30*b*Power(1/rhoa,0.6666666666666666),2)*
|
---|
3371 | Power(rhoa,4.333333333333333))) +
|
---|
3372 | (Power(beta,2)*Log(1 + (2*alpha*
|
---|
3373 | (-(alpha*e4ec*Power(mdr,4)*Power(Pi/6.,0.6666666666666666))/
|
---|
3374 | (32.*beta*(-1 + e4ec)*Power(rhoa,4.666666666666667)) +
|
---|
3375 | (Power(mdr,2)*Power(Pi/6.,0.3333333333333333))/
|
---|
3376 | (8.*Power(rhoa,2.3333333333333335))))/
|
---|
3377 | (beta + (Power(alpha,2)*e8ec*Power(mdr,4)*
|
---|
3378 | Power(Pi/6.,0.6666666666666666))/
|
---|
3379 | (16.*beta*Power(-1 + e4ec,2)*Power(rhoa,4.666666666666667)) -
|
---|
3380 | (alpha*e4ec*Power(mdr,2)*Power(Pi/6.,0.3333333333333333))/
|
---|
3381 | (4.*(-1 + e4ec)*Power(rhoa,2.3333333333333335)))))/(4.*alpha)
|
---|
3382 | ;
|
---|
3383 |
|
---|
3384 | return result;
|
---|
3385 | }
|
---|
3386 |
|
---|
3387 | void
|
---|
3388 | PW91CFunctional::point(const PointInputData &id, PointOutputData &od)
|
---|
3389 | {
|
---|
3390 | double ec_local, dec_local_rs, dec_local_zeta;
|
---|
3391 | local_->point_lc(id, od, ec_local, dec_local_rs, dec_local_zeta);
|
---|
3392 | double rho = id.a.rho+id.b.rho;
|
---|
3393 | double rho_13 = pow(rho, 1./3.);
|
---|
3394 | double rho_43 = rho*rho_13;
|
---|
3395 | double rs = 0.62035049089940009/rho_13;
|
---|
3396 | double z = (id.a.rho - id.b.rho)/rho;
|
---|
3397 | double gamma = sqrt(id.a.gamma + id.b.gamma + 2*id.gamma_ab);
|
---|
3398 |
|
---|
3399 | double pwc, dpwc_drs, dpwc_dg;
|
---|
3400 |
|
---|
3401 | if (rho < MIN_DENSITY) return;
|
---|
3402 | if (gamma < MIN_SQRTGAMMA) {
|
---|
3403 | if (!compute_potential_) return;
|
---|
3404 |
|
---|
3405 | double limit_dpwc_dgaa, limit_dpwc_dgab;
|
---|
3406 |
|
---|
3407 | double ga0;
|
---|
3408 | ga0 = 6.2337093539550051e7*b*c_x*pow(rs,7.0)+(6233709.3539550053
|
---|
3409 | *c_x*d-(4363.596547768504*b))*pow(rs,6.0)+(6233709.3539550053
|
---|
3410 | *c*c_x-(4363.596547768504*a))*pow(rs,5.0)+(6233709.3539550053
|
---|
3411 | *c_x-11205.715934669519)*pow(rs,4.0);
|
---|
3412 | limit_dpwc_dgaa = -((1.0*(1.4645918875615231*ga0*pow(z+1.0,
|
---|
3413 | 0.66666666666666663)+1.4645918875615231*ga0*pow(-(1.0*z)+
|
---|
3414 | 1.0,0.66666666666666663)))/pow(1.8929525610284735e7*b*pow(rs,
|
---|
3415 | 3.0)+1892952.5610284733*d*pow(rs,2.0)+1892952.5610284733*c*
|
---|
3416 | rs+1892952.5610284733,1.0));
|
---|
3417 |
|
---|
3418 | double gab0;
|
---|
3419 | gab0 = 6.2337093539550051e7*b*c_x*pow(rs,7.0)+(
|
---|
3420 | 6233709.3539550053*c_x*d-(4363.596547768504*b))*pow(rs,6.0)+(
|
---|
3421 | 6233709.3539550053*c*c_x-(4363.596547768504*a))*pow(rs,5.0)+(
|
---|
3422 | 6233709.3539550053*c_x-11205.715934669519)*pow(rs,4.0);
|
---|
3423 | limit_dpwc_dgab = -((1.0*(1.4645918875615231*gab0*pow(z+1.0,
|
---|
3424 | 0.66666666666666663)+1.4645918875615231*gab0*pow(-(1.0*z)+
|
---|
3425 | 1.0,0.66666666666666663)))/pow(9464762.8051423673*b*pow(rs,
|
---|
3426 | 3.0)+946476.28051423666*d*pow(rs,2.0)+946476.28051423666*c*
|
---|
3427 | rs+946476.28051423666,1.0));
|
---|
3428 |
|
---|
3429 | od.df_dgamma_aa += limit_dpwc_dgaa;
|
---|
3430 | od.df_dgamma_bb += limit_dpwc_dgaa;
|
---|
3431 | od.df_dgamma_ab += limit_dpwc_dgab;
|
---|
3432 |
|
---|
3433 | return;
|
---|
3434 | }
|
---|
3435 |
|
---|
3436 | double e0,e1,e2,e3,e4,e5,e6,e7,e8,e9,e10,e11,e12,e13,e14,e15;
|
---|
3437 | e2 = rs*rs; // pow(rs,2.0);
|
---|
3438 | e0 = e2*rs; // pow(rs,3.0);
|
---|
3439 | e1 = e0*e0*rs; // pow(rs,7.0);
|
---|
3440 | e3 = pow(z+1.0,2./3.)+pow(-z+1.0, 2./3.);
|
---|
3441 | e4 = gamma*gamma; // pow(gamma,2.0);
|
---|
3442 | e5 = e3*e3; // pow(e3,2.0);
|
---|
3443 | e6 = c_c0*c_c0; // pow(c_c0,2.0);
|
---|
3444 | e7 = e3*e3*e3; // pow(e3,3.0);
|
---|
3445 | e8 = 1./c_c0; // 1/pow(c_c0,1.0);
|
---|
3446 | e9 = 1/e5;
|
---|
3447 | e10 = 1/e6;
|
---|
3448 | e11 = exp(-0.064451410964169495*alpha*e10*ec_local/e7)-1.0;
|
---|
3449 | e12 = 1./e11; // 1/pow(e11,1.0);
|
---|
3450 | e13 = e1*e1; // pow(rs,14.0);
|
---|
3451 | e14 = 1/(e7*e3); // pow(e3,4.0);
|
---|
3452 | e15 = e4*e4; // pow(gamma,4.0);
|
---|
3453 | pwc = (0.238732414637843*((15.515564128703568*e6*e7*log((
|
---|
3454 | 0.12693641219540738*alpha*e8*(6.5448368524534253*alpha*e8*e13*
|
---|
3455 | e14*e12*e15+7.18052672676657*e1*e9*e4))/((
|
---|
3456 | 0.83077810845472067*alpha*alpha*e10*e13*e14*e15)/(e11*e11)
|
---|
3457 | +0.91147030036898102*alpha*e8*e1*e9*e12*e4+1.0)+
|
---|
3458 | 1.0))/alpha+(14.141975896783627*e1*((0.001*(b*e2+a
|
---|
3459 | *rs+2.568))/(10.0*b*e0+d*e2+c*rs+1.0)-(
|
---|
3460 | 1.4285714285714286*c_x)-(1.0*c_c0))*e3*e4)*exp(-
|
---|
3461 | 29.773894288156065*e1*rs*e5*e4)))/e0;
|
---|
3462 |
|
---|
3463 | if (compute_potential_) {
|
---|
3464 |
|
---|
3465 | double drs_drhoa = -0.20678349696646667/rho_43;
|
---|
3466 | double drs_drhob = drs_drhoa;
|
---|
3467 |
|
---|
3468 | double r0,r1,r2,r3,r4,r5,r6,r7,r8,r9,r10,r11,r12,r13,r14,r15,r16,r17,r18
|
---|
3469 | ,r19,r20,r21,r22,r23,r24,r25,r26,r27,r28,r29,r30,r31,r32,r33,t0;
|
---|
3470 |
|
---|
3471 | r0 = pow(rs,3.0);
|
---|
3472 | r1 = 1/pow(alpha,1.0);
|
---|
3473 | r2 = pow(c_c0,2.0);
|
---|
3474 | r3 = pow(z+1.0,0.66666666666666663)+pow(-(1.0*z)+1.0,
|
---|
3475 | 0.66666666666666663);
|
---|
3476 | r4 = pow(r3,3.0);
|
---|
3477 | r5 = 1/pow(c_c0,1.0);
|
---|
3478 | r6 = pow(rs,7.0);
|
---|
3479 | r7 = pow(r3,2.0);
|
---|
3480 | r8 = 1/r7;
|
---|
3481 | r9 = 1/r2;
|
---|
3482 | r10 = exp(-0.064451410964169495*alpha*r9*ec_local*1.0/r4);
|
---|
3483 | r11 = r10-1.0;
|
---|
3484 | r12 = 1/pow(r11,1.0);
|
---|
3485 | r13 = pow(gamma,2.0);
|
---|
3486 | r14 = pow(alpha,2.0);
|
---|
3487 | r15 = pow(rs,14.0);
|
---|
3488 | r16 = 1/pow(r3,4.0);
|
---|
3489 | r17 = 1/pow(r11,2.0);
|
---|
3490 | r18 = pow(gamma,4.0);
|
---|
3491 | r19 = 0.83077810845472067*r14*r9*r15*r16*r17*r18+
|
---|
3492 | 0.91147030036898102*alpha*r5*r6*r8*r12*r13+1.0;
|
---|
3493 | r20 = 1/pow(r19,1.0);
|
---|
3494 | r21 = 6.5448368524534253*alpha*r5*r15*r16*r12*r18+
|
---|
3495 | 7.18052672676657*r6*r8*r13;
|
---|
3496 | r22 = 0.12693641219540738*alpha*r5*r20*r21+1.0;
|
---|
3497 | r23 = pow(rs,6.0);
|
---|
3498 | r24 = 1/pow(c_c0,3.0);
|
---|
3499 | r25 = dec_local_rs;
|
---|
3500 | r26 = pow(rs,13.0);
|
---|
3501 | r27 = 1/pow(r3,7.0);
|
---|
3502 | r28 = pow(rs,2.0);
|
---|
3503 | r29 = b*r28+a*rs+2.568;
|
---|
3504 | r30 = 10.0*b*r0+d*r28+c*rs+1.0;
|
---|
3505 | r31 = 1/pow(r30,1.0);
|
---|
3506 | r32 = exp(-29.773894288156065*pow(rs,8.0)*r7*r13);
|
---|
3507 | r33 = 0.001*r29*r31-(1.4285714285714286*c_x)-(1.0*c_c0);
|
---|
3508 | t0 = -((0.71619724391352901*(15.515564128703568*r1*r2*r4*log
|
---|
3509 | (r22)+14.141975896783627*r6*r33*r3*r13*r32))/pow(rs,4.0));
|
---|
3510 | dpwc_drs = t0+(0.238732414637843*(-(3368.493563011893*r15*
|
---|
3511 | r33*r4*r18*r32)+98.993831277485398*r23*r33*r3*r13*r32+
|
---|
3512 | 14.141975896783627*r6*(0.001*(2.0*b*rs+a)*r31-((0.001*
|
---|
3513 | r29*(30.0*b*r28+2.0*d*rs+c))/pow(r30,2.0)))*r3*r13*r32+(
|
---|
3514 | 15.515564128703568*r1*r2*r4*(0.12693641219540738*alpha*r5*
|
---|
3515 | r20*(0.42182396967091729*r14*r24*r25*r15*r27*r17*r10*r18+
|
---|
3516 | 91.627715934347961*alpha*r5*r26*r16*r12*r18+
|
---|
3517 | 50.263687087365994*r23*r8*r13)-((0.12693641219540738*alpha*
|
---|
3518 | r5*r21*((0.10708964257610123*pow(alpha,3.0)*r25*r15*r27*r10
|
---|
3519 | *r18)/pow(c_c0,4.0)/pow(r11,3.0)+11.630893518366092*r14*
|
---|
3520 | r9*r26*r16*r17*r18+(0.058745546910716179*r14*r24*r25*r6*r17*
|
---|
3521 | r10*r13)/pow(r3,5.0)+6.3802921025828665*alpha*r5*r23*r8*r12
|
---|
3522 | *r13))/pow(r19,2.0))))/pow(r22,1.0)))/r0;
|
---|
3523 |
|
---|
3524 | if (id.a.rho > MIN_DENSITY && id.b.rho > MIN_DENSITY) {
|
---|
3525 | double dpwc_dz;
|
---|
3526 |
|
---|
3527 | double z0,z1,z2,z3,z4,z5,z6,z7,z8,z9,z10,z11,z12,z13,z14,z15,z16,z17
|
---|
3528 | ,z18 ,z19,z20,z21,z22,z23,z24,z25,z26,z27,z28,z29,z30,z31;
|
---|
3529 | z0 = pow(rs,3.0);
|
---|
3530 | z1 = 1/pow(alpha,1.0);
|
---|
3531 | z2 = pow(c_c0,2.0);
|
---|
3532 | z3 = -(1.0*z)+1.0;
|
---|
3533 | z4 = z+1.0;
|
---|
3534 | z5 = pow(z4,0.66666666666666663)+pow(z3,0.66666666666666663);
|
---|
3535 | z6 = pow(z5,3.0);
|
---|
3536 | z7 = 1/pow(c_c0,1.0);
|
---|
3537 | z8 = pow(rs,7.0);
|
---|
3538 | z9 = pow(z5,2.0);
|
---|
3539 | z10 = 1/z9;
|
---|
3540 | z11 = 1/z2;
|
---|
3541 | z12 = 1/z6;
|
---|
3542 | z13 = exp(-0.064451410964169495*alpha*z11*ec_local*z12);
|
---|
3543 | z14 = z13-1.0;
|
---|
3544 | z15 = 1/pow(z14,1.0);
|
---|
3545 | z16 = pow(gamma,2.0);
|
---|
3546 | z17 = pow(alpha,2.0);
|
---|
3547 | z18 = pow(rs,14.0);
|
---|
3548 | z19 = 1/pow(z5,4.0);
|
---|
3549 | z20 = 1/pow(z14,2.0);
|
---|
3550 | z21 = pow(gamma,4.0);
|
---|
3551 | z22 = 0.83077810845472067*z17*z11*z18*z19*z20*z21+
|
---|
3552 | 0.91147030036898102*alpha*z7*z8*z10*z15*z16+1.0;
|
---|
3553 | z23 = 1/pow(z22,1.0);
|
---|
3554 | z24 = 6.5448368524534253*alpha*z7*z18*z19*z15*z21+
|
---|
3555 | 7.18052672676657*z8*z10*z16;
|
---|
3556 | z25 = 0.12693641219540738*alpha*z7*z23*z24+1.0;
|
---|
3557 | z26 = 0.66666666666666663/pow(z4,0.33333333333333331)-(
|
---|
3558 | 0.66666666666666663/pow(z3,0.33333333333333331));
|
---|
3559 | z27 = -(0.064451410964169495*alpha*z11*dec_local_zeta*z12)
|
---|
3560 | +0.19335423289250847*alpha*z11*ec_local*z26*z19;
|
---|
3561 | z28 = 1/pow(z5,5.0);
|
---|
3562 | z29 = pow(rs,2.0);
|
---|
3563 | z30 = (0.001*(b*z29+a*rs+2.5680000000000001))/pow(10.0*b*z0+d*
|
---|
3564 | z29+c*rs+1.0,1.0)-(1.4285714285714286*c_x)-(1.0*c_c0);
|
---|
3565 | z31 = exp(-29.773894288156065*pow(rs,8.0)*z9*z16);
|
---|
3566 | dpwc_dz = (0.238732414637843*(46.546692386110699*z1*z2*z26*z9*
|
---|
3567 | log(z25)-(842.12339075297325*pow(rs,15.0)*z30*z26*z9*z21*z31)+
|
---|
3568 | 14.141975896783627*z8*z30*z26*z16*z31+(15.515564128703568*z1*z2
|
---|
3569 | *z6*(0.12693641219540738*alpha*z7*z23*(-(6.5448368524534253*
|
---|
3570 | alpha*z7*z18*z19*z27*z20*z13*z21)-(26.179347409813701*alpha*z7*
|
---|
3571 | z18*z26*z28*z15*z21)-(14.36105345353314*z8*z26*z12*z16))-((
|
---|
3572 | 0.12693641219540738*alpha*z7*z24*(-((1.6615562169094413*z17*z11
|
---|
3573 | *z18*z19*z27*z13*z21)/pow(z14,3.0))-(3.3231124338188827*z17*z11
|
---|
3574 | *z18*z26*z28*z20*z21)-(0.91147030036898102*alpha*z7*z8*z10*z27*
|
---|
3575 | z20*z13*z16)-(1.822940600737962*alpha*z7*z8*z26*z12*z15*z16)))/
|
---|
3576 | pow(z22,2.0))))/pow(z25,1.0)))/z0;
|
---|
3577 |
|
---|
3578 | double dz_drhoa = ( 1. - (id.a.rho-id.b.rho)/rho)/rho;
|
---|
3579 | double dz_drhob = (-1. - (id.a.rho-id.b.rho)/rho)/rho;
|
---|
3580 |
|
---|
3581 | od.df_drho_a += dpwc_drs * drs_drhoa + dpwc_dz * dz_drhoa;
|
---|
3582 | od.df_drho_b += dpwc_drs * drs_drhob + dpwc_dz * dz_drhob;
|
---|
3583 | }
|
---|
3584 | else if (id.a.rho > MIN_DENSITY) {
|
---|
3585 | // df_drho_b diverges
|
---|
3586 | double dzeta_drhoa = 1./rho * (1. - z);
|
---|
3587 | double drs_drhoa = -rs/(3.*rho);
|
---|
3588 | double dec_local_drhoa = dec_local_rs*drs_drhoa
|
---|
3589 | + dec_local_zeta*dzeta_drhoa;
|
---|
3590 | od.df_drho_a += limit_df_drhoa(id.a.rho, gamma,
|
---|
3591 | ec_local, dec_local_drhoa);
|
---|
3592 | }
|
---|
3593 | else if (id.b.rho > MIN_DENSITY) {
|
---|
3594 | // df_drho_a diverges
|
---|
3595 | double dzeta_drhob = 1./rho * (-1. - z);
|
---|
3596 | double drs_drhob = -rs/(3.*rho);
|
---|
3597 | double dec_local_drhob = dec_local_rs*drs_drhob
|
---|
3598 | + dec_local_zeta*dzeta_drhob;
|
---|
3599 | od.df_drho_b += limit_df_drhoa(id.b.rho, gamma,
|
---|
3600 | ec_local, dec_local_drhob);
|
---|
3601 | }
|
---|
3602 |
|
---|
3603 | double g0,g1,g2,g3,g4,g5,g6,g7,g8,g9,g10,g11,g12,g13,g14,g15,g16,g17,g18
|
---|
3604 | ,g19,g20,g21,g22,g23;
|
---|
3605 | g0 = pow(rs,3.0);
|
---|
3606 | g1 = pow(c_c0,2.0);
|
---|
3607 | g2 = pow(z+1.0,0.66666666666666663)+pow(-(1.0*z)+1.0,
|
---|
3608 | 0.66666666666666663);
|
---|
3609 | g3 = pow(g2,3.0);
|
---|
3610 | g4 = 1/pow(c_c0,1.0);
|
---|
3611 | g5 = pow(rs,7.0);
|
---|
3612 | g6 = pow(g2,2.0);
|
---|
3613 | g7 = 1/g6;
|
---|
3614 | g8 = pow(rs,14.0);
|
---|
3615 | g9 = 1/pow(g2,4.0);
|
---|
3616 | g10 = 1/g1;
|
---|
3617 | g11 = exp(-0.064451410964169495*alpha*g10*ec_local*1.0/g3)-
|
---|
3618 | 1.0;
|
---|
3619 | g12 = 1/pow(g11,1.0);
|
---|
3620 | g13 = pow(gamma,3.0);
|
---|
3621 | g14 = pow(gamma,2.0);
|
---|
3622 | g15 = pow(alpha,2.0);
|
---|
3623 | g16 = 1/pow(g11,2.0);
|
---|
3624 | g17 = pow(gamma,4.0);
|
---|
3625 | g18 = 0.83077810845472067*g15*g10*g8*g9*g16*g17+
|
---|
3626 | 0.91147030036898102*alpha*g4*g5*g7*g12*g14+1.0;
|
---|
3627 | g19 = 1/pow(g18,1.0);
|
---|
3628 | g20 = 6.5448368524534253*alpha*g4*g8*g9*g12*g17+7.18052672676657
|
---|
3629 | *g5*g7*g14;
|
---|
3630 | g21 = pow(rs,2.0);
|
---|
3631 | g22 = (0.001*(b*g21+a*rs+2.5680000000000001))/pow(10.0*b*g0+d*
|
---|
3632 | g21+c*rs+1.0,1.0)-(1.4285714285714286*c_x)-(1.0*c_c0);
|
---|
3633 | g23 = exp(-29.773894288156065*pow(rs,8.0)*g6*g14);
|
---|
3634 | dpwc_dg = (0.238732414637843*(-(842.12339075297325*pow(rs,15.0
|
---|
3635 | )*g22*g3*g13*g23)+28.283951793567255*g5*g22*g2*gamma*g23+(
|
---|
3636 | 15.515564128703568*g1*g3*(-((0.12693641219540738*alpha*g4*(
|
---|
3637 | 3.3231124338188827*g15*g10*g8*g9*g16*g13+1.822940600737962*
|
---|
3638 | alpha*g4*g5*g7*g12*gamma)*g20)/pow(g18,2.0))+
|
---|
3639 | 0.12693641219540738*alpha*g4*(26.179347409813701*alpha*g4*g8*g9
|
---|
3640 | *g12*g13+14.36105345353314*g5*g7*gamma)*g19))/pow(alpha,1.0)/
|
---|
3641 | pow(0.12693641219540738*alpha*g4*g19*g20+1.0,1.0)))/g0;
|
---|
3642 |
|
---|
3643 | od.df_dgamma_aa += dpwc_dg / (2 * gamma);
|
---|
3644 | od.df_dgamma_bb += dpwc_dg / (2 * gamma);
|
---|
3645 | od.df_dgamma_ab += dpwc_dg / gamma;
|
---|
3646 | }
|
---|
3647 |
|
---|
3648 | od.energy += pwc;
|
---|
3649 | }
|
---|
3650 |
|
---|
3651 | /////////////////////////////////////////////////////////////////////////////
|
---|
3652 | // PW91XFunctional
|
---|
3653 |
|
---|
3654 | static ClassDesc PW91XFunctional_cd(
|
---|
3655 | typeid(PW91XFunctional),"PW91XFunctional",1,"public DenFunctional",
|
---|
3656 | 0, create<PW91XFunctional>, create<PW91XFunctional>);
|
---|
3657 |
|
---|
3658 | PW91XFunctional::PW91XFunctional(StateIn& s):
|
---|
3659 | SavableState(s),
|
---|
3660 | DenFunctional(s)
|
---|
3661 | {
|
---|
3662 | s.get(a);
|
---|
3663 | s.get(b);
|
---|
3664 | s.get(c);
|
---|
3665 | s.get(d);
|
---|
3666 | s.get(a_x);
|
---|
3667 | }
|
---|
3668 |
|
---|
3669 | PW91XFunctional::PW91XFunctional()
|
---|
3670 | {
|
---|
3671 | init_constants();
|
---|
3672 | }
|
---|
3673 |
|
---|
3674 | PW91XFunctional::PW91XFunctional(const Ref<KeyVal>& keyval):
|
---|
3675 | DenFunctional(keyval)
|
---|
3676 | {
|
---|
3677 | init_constants();
|
---|
3678 | a_x = keyval->doublevalue("a_x", KeyValValuedouble(a_x));
|
---|
3679 | a = keyval->doublevalue("a", KeyValValuedouble(a));
|
---|
3680 | b = keyval->doublevalue("b", KeyValValuedouble(b));
|
---|
3681 | c = keyval->doublevalue("c", KeyValValuedouble(c));
|
---|
3682 | d = keyval->doublevalue("d", KeyValValuedouble(d));
|
---|
3683 | }
|
---|
3684 |
|
---|
3685 | PW91XFunctional::~PW91XFunctional()
|
---|
3686 | {
|
---|
3687 | }
|
---|
3688 |
|
---|
3689 | void
|
---|
3690 | PW91XFunctional::init_constants()
|
---|
3691 | {
|
---|
3692 | a = 0.19645;
|
---|
3693 | b = 7.7956;
|
---|
3694 | c = 0.2743;
|
---|
3695 | // the PW91 paper had d = 0.1508
|
---|
3696 | // PBE.f has the following
|
---|
3697 | d = 0.15084;
|
---|
3698 | // a_x -(3/4)*(3/pi)^(1/3)
|
---|
3699 | // the following rounded a_x appears in PBE.f
|
---|
3700 | a_x = -0.7385588;
|
---|
3701 | }
|
---|
3702 |
|
---|
3703 | void
|
---|
3704 | PW91XFunctional::save_data_state(StateOut& s)
|
---|
3705 | {
|
---|
3706 | DenFunctional::save_data_state(s);
|
---|
3707 | s.put(a);
|
---|
3708 | s.put(b);
|
---|
3709 | s.put(c);
|
---|
3710 | s.put(d);
|
---|
3711 | s.put(a_x);
|
---|
3712 | }
|
---|
3713 |
|
---|
3714 | int
|
---|
3715 | PW91XFunctional::need_density_gradient()
|
---|
3716 | {
|
---|
3717 | return 1;
|
---|
3718 | }
|
---|
3719 |
|
---|
3720 | void
|
---|
3721 | PW91XFunctional::spin_contrib(const PointInputData::SpinData &i,
|
---|
3722 | double &pw, double &dpw_dr, double &dpw_dg)
|
---|
3723 | {
|
---|
3724 | double rho = i.rho;
|
---|
3725 | double gamma = i.gamma;
|
---|
3726 | const double pi = M_PI;
|
---|
3727 |
|
---|
3728 | if (rho < MIN_DENSITY) {
|
---|
3729 | pw = 0.;
|
---|
3730 | dpw_dr = 0.;
|
---|
3731 | dpw_dg = 0.; // really -inf
|
---|
3732 | return;
|
---|
3733 | }
|
---|
3734 | if (gamma < MIN_GAMMA) {
|
---|
3735 | double rho_43 = rho * i.rho_13; // rho^(4/3)
|
---|
3736 | // 2^(1/3) * a_x * rho^(4/3)
|
---|
3737 | pw = 1.2599210498948732 * a_x * rho_43;
|
---|
3738 | // (4/3) 2^(1/3) a_x rho^(1/3)
|
---|
3739 | dpw_dr = 1.6798947331931642 * a_x * i.rho_13;
|
---|
3740 | // (6^(2/3) a_x / (24 3^(1/3) pi^4/3)) (d - c) / rho^(4/3)
|
---|
3741 | dpw_dg = - 0.020732388737701564 * a_x * (d - c) / rho_43;
|
---|
3742 | return;
|
---|
3743 | }
|
---|
3744 |
|
---|
3745 | // this has been generated by macsyma and modified
|
---|
3746 | double t0,t1,t2,t3,t4,t5,t6,t7,t8,t9;
|
---|
3747 | t0 = 1.8171205928321397; // pow(6,1.0/3.0);
|
---|
3748 | t1 = rho * i.rho_13; // pow(rho,4.0/3.0);
|
---|
3749 | t2 = 1/t0;
|
---|
3750 | t3 = 0.46619407703541166; // 1/pow(pi,2.0/3.0);
|
---|
3751 | t4 = 1/t1;
|
---|
3752 | t5 = sqrt(gamma);
|
---|
3753 | t6 = (t2*a*t3*t4*asinh((t2*b*t3*t4*t5)/2.0)*t5)/2.0;
|
---|
3754 | t7 = 0.30285343213868998; // 1/pow(6,2.0/3.0);
|
---|
3755 | t8 = 0.21733691746289932; // 1/pow(pi,4.0/3.0);
|
---|
3756 | t9 = t4*t4; // 1/pow(rho,8.0/3.0);
|
---|
3757 | pw = (t0*a_x*t1*((t7*t8*t9*gamma*(-(d*exp(-25.0*t7*t8*t9*gamma))+c)
|
---|
3758 | )/4.0+t6+1))/pow(3,1.0/3.0)/((t2*pow(gamma,2.0))/24000.0/pow(pi,8.0
|
---|
3759 | /3.0)/pow(rho,16.0/3.0)+t6+1);
|
---|
3760 | if (compute_potential_) {
|
---|
3761 | double r0,r1,r2,r3,r4,r5,r6,r7,r8,r9,r10,r11,r12,r13,r14,r15,r16,r17,r18
|
---|
3762 | ,r19,r20,r21,r22;
|
---|
3763 | r0 = 0.69336127435063466; // 1/pow(3,1.0/3.0);
|
---|
3764 | r1 = t0; // pow(6,1.0/3.0);
|
---|
3765 | r2 = t1; // pow(rho,4.0/3.0);
|
---|
3766 | r3 = 1/r1;
|
---|
3767 | r4 = t3; // 1/pow(pi,2.0/3.0);
|
---|
3768 | r5 = 1/r2;
|
---|
3769 | r6 = t5; // sqrt(gamma);
|
---|
3770 | r7 = asinh((r3*b*r4*r5*r6)/2.0);
|
---|
3771 | r8 = (r3*a*r4*r5*r7*r6)/2.0;
|
---|
3772 | r9 = 1/pow(pi,8.0/3.0);
|
---|
3773 | r10 = gamma*gamma; // pow(gamma,2.0);
|
---|
3774 | r11 = (r3*r9*r10)/24000.0/pow(rho,16.0/3.0)+r8+1;
|
---|
3775 | r12 = 1/r11;
|
---|
3776 | r13 = -((2.0/3.0*r3*a*r4*r7*r6)/pow(rho,7.0/3.0));
|
---|
3777 | r14 = 1/pow(6,2.0/3.0);
|
---|
3778 | r15 = 1/pow(pi,4.0/3.0);
|
---|
3779 | r16 = 1/pow(rho,11.0/3.0);
|
---|
3780 | r17 = 1/pow(rho,8.0/3.0);
|
---|
3781 | r18 = -((1.0/3.0*r14*a*b*r15*r16*gamma)/sqrt((r14*pow(b,2.0)*r15*r17
|
---|
3782 | *gamma)/4.0+1));
|
---|
3783 | r19 = 1/pow(rho,19.0/3.0);
|
---|
3784 | r20 = exp(-25.0*r14*r15*r17*gamma);
|
---|
3785 | r21 = -(d*r20)+c;
|
---|
3786 | r22 = (r14*r15*r17*gamma*r21)/4.0+r8+1;
|
---|
3787 | dpw_dr = -((r0*r1*a_x*r2*(r18-(1.0/4500.0*r3*r9*r19*r10)+r13)*r22)/
|
---|
3788 | pow(r11,2.0))+4.0/3.0*r0*r1*a_x*pow(rho,1.0/3.0)*r12*r22+r0*r1*a_x*
|
---|
3789 | r2*r12*(-(2.0/3.0*r14*r15*r16*gamma*r21)-(25.0/9.0*r3*d*r9*r19*r10*
|
---|
3790 | r20)+r18+r13);
|
---|
3791 | double g0,g1,g2,g3,g4,g5,g6,g7,g8,g9,g10,g11,g12,g13,g14,g15,g16,g17,g18;
|
---|
3792 | g0 = 1/pow(3,1.0/3.0);
|
---|
3793 | g1 = pow(6,1.0/3.0);
|
---|
3794 | g2 = pow(rho,4.0/3.0);
|
---|
3795 | g3 = 1/g1;
|
---|
3796 | g4 = 1/pow(pi,2.0/3.0);
|
---|
3797 | g5 = 1/g2;
|
---|
3798 | g6 = sqrt(gamma);
|
---|
3799 | g7 = asinh((g3*b*g4*g5*g6)/2.0);
|
---|
3800 | g8 = (g3*a*g4*g5*g7*g6)/2.0;
|
---|
3801 | g9 = 1/pow(pi,8.0/3.0);
|
---|
3802 | g10 = 1/pow(rho,16.0/3.0);
|
---|
3803 | g11 = (g3*g9*g10*pow(gamma,2.0))/24000.0+g8+1;
|
---|
3804 | g12 = (g3*a*g4*g5*g7)/4.0/g6;
|
---|
3805 | g13 = 1/pow(6,2.0/3.0);
|
---|
3806 | g14 = 1/pow(pi,4.0/3.0);
|
---|
3807 | g15 = 1/pow(rho,8.0/3.0);
|
---|
3808 | g16 = (g13*a*b*g14*g15)/8.0/sqrt((g13*pow(b,2.0)*g14*g15*gamma)/4.0+
|
---|
3809 | 1);
|
---|
3810 | g17 = exp(-25.0*g13*g14*g15*gamma);
|
---|
3811 | g18 = -(d*g17)+c;
|
---|
3812 | dpw_dg = -((g0*g1*a_x*g2*(g16+(g3*g9*g10*gamma)/12000.0+g12)*((g13*
|
---|
3813 | g14*g15*gamma*g18)/4.0+g8+1))/pow(g11,2.0))+(g0*g1*a_x*g2*((g13*g14
|
---|
3814 | *g15*g18)/4.0+25.0/24.0*g3*d*g9*g10*gamma*g17+g16+g12))/g11;
|
---|
3815 | }
|
---|
3816 | }
|
---|
3817 |
|
---|
3818 | void
|
---|
3819 | PW91XFunctional::point(const PointInputData &id,
|
---|
3820 | PointOutputData &od)
|
---|
3821 | {
|
---|
3822 | od.zero();
|
---|
3823 |
|
---|
3824 | double mpw, dmpw_dr, dmpw_dg;
|
---|
3825 |
|
---|
3826 | // alpha
|
---|
3827 | spin_contrib(id.a, mpw, dmpw_dr, dmpw_dg);
|
---|
3828 | od.energy = mpw;
|
---|
3829 | od.df_drho_a = dmpw_dr;
|
---|
3830 | od.df_dgamma_aa = dmpw_dg;
|
---|
3831 |
|
---|
3832 | // beta
|
---|
3833 | if (spin_polarized_) {
|
---|
3834 | spin_contrib(id.b, mpw, dmpw_dr, dmpw_dg);
|
---|
3835 | od.energy += mpw;
|
---|
3836 | od.df_drho_b = dmpw_dr;
|
---|
3837 | od.df_dgamma_bb = dmpw_dg;
|
---|
3838 | }
|
---|
3839 | else {
|
---|
3840 | od.energy += mpw;
|
---|
3841 | od.df_drho_b = od.df_drho_a;
|
---|
3842 | od.df_dgamma_bb = od.df_dgamma_aa;
|
---|
3843 | }
|
---|
3844 | }
|
---|
3845 |
|
---|
3846 | /////////////////////////////////////////////////////////////////////////////
|
---|
3847 | // PBEXFunctional
|
---|
3848 |
|
---|
3849 | static ClassDesc PBEXFunctional_cd(
|
---|
3850 | typeid(PBEXFunctional),"PBEXFunctional",1,"public DenFunctional",
|
---|
3851 | 0, create<PBEXFunctional>, create<PBEXFunctional>);
|
---|
3852 |
|
---|
3853 | PBEXFunctional::PBEXFunctional(StateIn& s):
|
---|
3854 | SavableState(s),
|
---|
3855 | DenFunctional(s)
|
---|
3856 | {
|
---|
3857 | s.get(mu);
|
---|
3858 | s.get(kappa);
|
---|
3859 | }
|
---|
3860 |
|
---|
3861 | PBEXFunctional::PBEXFunctional()
|
---|
3862 | {
|
---|
3863 | init_constants();
|
---|
3864 | }
|
---|
3865 |
|
---|
3866 | PBEXFunctional::PBEXFunctional(const Ref<KeyVal>& keyval):
|
---|
3867 | DenFunctional(keyval)
|
---|
3868 | {
|
---|
3869 | init_constants();
|
---|
3870 | mu = keyval->doublevalue("mu", KeyValValuedouble(mu));
|
---|
3871 | kappa = keyval->doublevalue("kappa", KeyValValuedouble(kappa));
|
---|
3872 | if (keyval->booleanvalue("revPBE")) {
|
---|
3873 | kappa = 1.245;
|
---|
3874 | }
|
---|
3875 | }
|
---|
3876 |
|
---|
3877 | PBEXFunctional::~PBEXFunctional()
|
---|
3878 | {
|
---|
3879 | }
|
---|
3880 |
|
---|
3881 | void
|
---|
3882 | PBEXFunctional::init_constants()
|
---|
3883 | {
|
---|
3884 | // in paper:
|
---|
3885 | // mu = 0.21951;
|
---|
3886 | // in PBE.F
|
---|
3887 | mu = 0.2195149727645171;
|
---|
3888 | kappa = 0.804;
|
---|
3889 | }
|
---|
3890 |
|
---|
3891 | void
|
---|
3892 | PBEXFunctional::save_data_state(StateOut& s)
|
---|
3893 | {
|
---|
3894 | DenFunctional::save_data_state(s);
|
---|
3895 | s.put(mu);
|
---|
3896 | s.put(kappa);
|
---|
3897 | }
|
---|
3898 |
|
---|
3899 | int
|
---|
3900 | PBEXFunctional::need_density_gradient()
|
---|
3901 | {
|
---|
3902 | return 1;
|
---|
3903 | }
|
---|
3904 |
|
---|
3905 | void
|
---|
3906 | PBEXFunctional::spin_contrib(const PointInputData::SpinData &i,
|
---|
3907 | double &pbex,
|
---|
3908 | double &dpbex_drhoa, double &dpbex_dgaa)
|
---|
3909 | {
|
---|
3910 | double rhoa = i.rho;
|
---|
3911 | double rhoa_13 = i.rho_13;
|
---|
3912 | double rhoa_43 = rhoa*rhoa_13;
|
---|
3913 | double gaa = i.gamma;
|
---|
3914 |
|
---|
3915 | if (rhoa < MIN_DENSITY) {
|
---|
3916 | pbex = 0;
|
---|
3917 | if (compute_potential_) {
|
---|
3918 | dpbex_drhoa = 0.0;
|
---|
3919 | dpbex_dgaa = 0.0; // really -inf
|
---|
3920 | }
|
---|
3921 | return;
|
---|
3922 | }
|
---|
3923 |
|
---|
3924 | if (gaa < MIN_GAMMA) {
|
---|
3925 | pbex = - 0.93052573634910019 * rhoa_43;
|
---|
3926 | if (compute_potential_) {
|
---|
3927 | dpbex_drhoa = -(1.2407009817988002*rhoa_13);
|
---|
3928 | dpbex_dgaa = -((0.015312087450269402*mu)/rhoa_43);
|
---|
3929 | }
|
---|
3930 | return;
|
---|
3931 | }
|
---|
3932 |
|
---|
3933 | double rhoa_83 = rhoa_43*rhoa_43;
|
---|
3934 |
|
---|
3935 | pbex = -(0.93052573634910007*(-(kappa/((
|
---|
3936 | 0.016455307846020562*gaa*mu)/kappa/rhoa_83+1.0))+kappa+1.0)*rhoa_43);
|
---|
3937 |
|
---|
3938 | if (compute_potential_) {
|
---|
3939 | double rhoa_73 = rhoa_43*rhoa;
|
---|
3940 | double r0;
|
---|
3941 | r0 = (0.016455307846020562*gaa*mu)/kappa/rhoa_83+1.0;
|
---|
3942 | dpbex_drhoa = -(1.2407009817988002*(-(kappa/r0)
|
---|
3943 | +kappa+1.0)*rhoa_13)+(0.040832233200718403*gaa*mu)/(r0*r0)/rhoa_73;
|
---|
3944 | dpbex_dgaa = -((0.015312087450269402*mu)/(r0*r0)/rhoa_43);
|
---|
3945 | }
|
---|
3946 | }
|
---|
3947 |
|
---|
3948 | void
|
---|
3949 | PBEXFunctional::point(const PointInputData &id,
|
---|
3950 | PointOutputData &od)
|
---|
3951 | {
|
---|
3952 | od.zero();
|
---|
3953 |
|
---|
3954 | double pbex, dpbex_dr, dpbex_dg;
|
---|
3955 |
|
---|
3956 | // alpha
|
---|
3957 | spin_contrib(id.a, pbex, dpbex_dr, dpbex_dg);
|
---|
3958 | od.energy = pbex;
|
---|
3959 | if (compute_potential_) {
|
---|
3960 | od.df_drho_a = dpbex_dr;
|
---|
3961 | od.df_dgamma_aa = dpbex_dg;
|
---|
3962 | }
|
---|
3963 |
|
---|
3964 | // beta
|
---|
3965 | if (spin_polarized_) {
|
---|
3966 | spin_contrib(id.b, pbex, dpbex_dr, dpbex_dg);
|
---|
3967 | od.energy += pbex;
|
---|
3968 | if (compute_potential_) {
|
---|
3969 | od.df_drho_b = dpbex_dr;
|
---|
3970 | od.df_dgamma_bb = dpbex_dg;
|
---|
3971 | }
|
---|
3972 | }
|
---|
3973 | else {
|
---|
3974 | od.energy += pbex;
|
---|
3975 | if (compute_potential_) {
|
---|
3976 | od.df_drho_b = od.df_drho_a;
|
---|
3977 | od.df_dgamma_bb = od.df_dgamma_aa;
|
---|
3978 | }
|
---|
3979 | }
|
---|
3980 | }
|
---|
3981 |
|
---|
3982 | /////////////////////////////////////////////////////////////////////////////
|
---|
3983 | // mPW91XFunctional
|
---|
3984 |
|
---|
3985 | static ClassDesc mPW91XFunctional_cd(
|
---|
3986 | typeid(mPW91XFunctional),"mPW91XFunctional",1,"public DenFunctional",
|
---|
3987 | 0, create<mPW91XFunctional>, create<mPW91XFunctional>);
|
---|
3988 |
|
---|
3989 | mPW91XFunctional::mPW91XFunctional(StateIn& s):
|
---|
3990 | SavableState(s),
|
---|
3991 | DenFunctional(s)
|
---|
3992 | {
|
---|
3993 | init_constants(mPW91);
|
---|
3994 | s.get(b);
|
---|
3995 | s.get(beta);
|
---|
3996 | s.get(c);
|
---|
3997 | s.get(d);
|
---|
3998 | s.get(x_d_coef);
|
---|
3999 | }
|
---|
4000 |
|
---|
4001 | mPW91XFunctional::mPW91XFunctional()
|
---|
4002 | {
|
---|
4003 | init_constants(mPW91);
|
---|
4004 | }
|
---|
4005 |
|
---|
4006 | mPW91XFunctional::mPW91XFunctional(mPW91XFunctional::Func f)
|
---|
4007 | {
|
---|
4008 | init_constants(f);
|
---|
4009 | }
|
---|
4010 |
|
---|
4011 | mPW91XFunctional::mPW91XFunctional(const Ref<KeyVal>& keyval):
|
---|
4012 | DenFunctional(keyval)
|
---|
4013 | {
|
---|
4014 | char *t = keyval->pcharvalue("constants");
|
---|
4015 | if (t) {
|
---|
4016 | if (!strcmp(t,"B88")) {
|
---|
4017 | init_constants(B88);
|
---|
4018 | }
|
---|
4019 | else if (!strcmp(t,"PW91")) {
|
---|
4020 | init_constants(PW91);
|
---|
4021 | }
|
---|
4022 | else if (!strcmp(t,"mPW91")) {
|
---|
4023 | init_constants(mPW91);
|
---|
4024 | }
|
---|
4025 | else {
|
---|
4026 | ExEnv::outn() << "mPW91XFunctional: bad \"constants\": " << t << endl;
|
---|
4027 | abort();
|
---|
4028 | }
|
---|
4029 | delete[] t;
|
---|
4030 | }
|
---|
4031 | else {
|
---|
4032 | init_constants(mPW91);
|
---|
4033 | }
|
---|
4034 | b = keyval->doublevalue("b", KeyValValuedouble(b));
|
---|
4035 | beta = keyval->doublevalue("beta", KeyValValuedouble(beta));
|
---|
4036 | c = keyval->doublevalue("c", KeyValValuedouble(c));
|
---|
4037 | d = keyval->doublevalue("d", KeyValValuedouble(d));
|
---|
4038 | x_d_coef = keyval->doublevalue("x_d_coef", KeyValValuedouble(x_d_coef));
|
---|
4039 | }
|
---|
4040 |
|
---|
4041 | mPW91XFunctional::~mPW91XFunctional()
|
---|
4042 | {
|
---|
4043 | }
|
---|
4044 |
|
---|
4045 | void
|
---|
4046 | mPW91XFunctional::init_constants(Func f)
|
---|
4047 | {
|
---|
4048 | a_x = -1.5*pow(3./(4.*M_PI), 1./3.);
|
---|
4049 | if (f == B88) {
|
---|
4050 | b = 0.0042;
|
---|
4051 | beta = 0.0042;
|
---|
4052 | c = 0.;
|
---|
4053 | d = 1.;
|
---|
4054 | x_d_coef = 0.;
|
---|
4055 | }
|
---|
4056 | else if (f == PW91) {
|
---|
4057 | b = 0.0042;
|
---|
4058 | beta = 5.*pow(36.*M_PI,-5./3.);
|
---|
4059 | c = 1.6455;
|
---|
4060 | d = 4.;
|
---|
4061 | x_d_coef = 1.e-6;
|
---|
4062 | }
|
---|
4063 | else {
|
---|
4064 | b = 0.0046;
|
---|
4065 | beta = 5.*pow(36.*M_PI,-5./3.);
|
---|
4066 | c = 1.6455;
|
---|
4067 | d = 3.73;
|
---|
4068 | x_d_coef = 1.e-6;
|
---|
4069 | }
|
---|
4070 | }
|
---|
4071 |
|
---|
4072 | void
|
---|
4073 | mPW91XFunctional::save_data_state(StateOut& s)
|
---|
4074 | {
|
---|
4075 | DenFunctional::save_data_state(s);
|
---|
4076 | s.put(b);
|
---|
4077 | s.put(beta);
|
---|
4078 | s.put(c);
|
---|
4079 | s.put(d);
|
---|
4080 | s.put(x_d_coef);
|
---|
4081 | }
|
---|
4082 |
|
---|
4083 | int
|
---|
4084 | mPW91XFunctional::need_density_gradient()
|
---|
4085 | {
|
---|
4086 | return 1;
|
---|
4087 | }
|
---|
4088 |
|
---|
4089 | void
|
---|
4090 | mPW91XFunctional::spin_contrib(const PointInputData::SpinData &i,
|
---|
4091 | double &mpw, double &dmpw_dr, double &dmpw_dg)
|
---|
4092 | {
|
---|
4093 | // this has been generated by macsyma
|
---|
4094 | double rho = i.rho;
|
---|
4095 | double gamma = i.gamma;
|
---|
4096 |
|
---|
4097 | if (rho < MIN_DENSITY) {
|
---|
4098 | mpw = 0.;
|
---|
4099 | dmpw_dr = 0.;
|
---|
4100 | dmpw_dg = 0.; // division by zero
|
---|
4101 | return;
|
---|
4102 | }
|
---|
4103 | if (gamma < MIN_GAMMA) {
|
---|
4104 | double rho_43 = rho * i.rho_13; // rho^(4/3)
|
---|
4105 | // 2^(1/3) * a_x * rho^(4/3)
|
---|
4106 | mpw = a_x * rho_43;
|
---|
4107 | // (4/3) a_x rho^(1/3)
|
---|
4108 | dmpw_dr = (4./3.) * a_x * i.rho_13;
|
---|
4109 | dmpw_dg = -beta/rho_43;
|
---|
4110 | return;
|
---|
4111 | }
|
---|
4112 |
|
---|
4113 | double t0,t1,t2,t3,t4,t5;
|
---|
4114 | t0 = pow(rho,4.0/3.0);
|
---|
4115 | t1 = 1/t0;
|
---|
4116 | t2 = sqrt(gamma);
|
---|
4117 | t3 = 1/pow(rho,4.0/3.0*d);
|
---|
4118 | t4 = pow(gamma,d/2.0);
|
---|
4119 | t5 = 1/pow(rho,8.0/3.0);
|
---|
4120 | mpw = t0*(-((-(((-beta+b)*t5*gamma)*exp(-c*t5*gamma))-(t3*x_d_coef*t4
|
---|
4121 | )+b*t5*gamma)/(-((t3*x_d_coef*t4)/a_x)+6*b*t1*asinh(t1*t2)*t2+1))+
|
---|
4122 | a_x);
|
---|
4123 | if (compute_potential_) {
|
---|
4124 | double r0,r1,r2,r3,r4,r5,r6,r7,r8,r9,r10,r11,r12,r13,r14;
|
---|
4125 | r0 = pow(rho,4.0/3.0);
|
---|
4126 | r1 = 1/r0;
|
---|
4127 | r2 = sqrt(gamma);
|
---|
4128 | r3 = asinh(r1*r2);
|
---|
4129 | r4 = 1/a_x;
|
---|
4130 | r5 = 1/pow(rho,4.0/3.0*d);
|
---|
4131 | r6 = pow(gamma,d/2.0);
|
---|
4132 | r7 = -(r4*r5*x_d_coef*r6)+6*b*r1*r3*r2+1;
|
---|
4133 | r8 = 1/r7;
|
---|
4134 | r9 = 1/pow(rho,8.0/3.0);
|
---|
4135 | r10 = -beta+b;
|
---|
4136 | r11 = exp(-c*r9*gamma);
|
---|
4137 | r12 = -(r10*r9*gamma*r11)-(r5*x_d_coef*r6)+b*r9*gamma;
|
---|
4138 | r13 = 1/pow(rho,11.0/3.0);
|
---|
4139 | r14 = pow(rho,-(4.0/3.0*d)-1);
|
---|
4140 | dmpw_dr = r0*(-(r8*(-((8.0/3.0*r10*c*pow(gamma,2.0)*r11)/pow(rho,
|
---|
4141 | 19.0/3.0))+8.0/3.0*r10*r13*gamma*r11+4.0/3.0*d*r14*x_d_coef*r6-(8.0
|
---|
4142 | /3.0*b*r13*gamma)))+((4.0/3.0*r4*d*r14*x_d_coef*r6-((8*b*r13*gamma)
|
---|
4143 | /sqrt(r9*gamma+1))-((8*b*r3*r2)/pow(rho,7.0/3.0)))*r12)/pow(r7,2.0)
|
---|
4144 | )+4.0/3.0*pow(rho,1.0/3.0)*(-(r8*r12)+a_x);
|
---|
4145 | double g0,g1,g2,g3,g4,g5,g6,g7,g8,g9,g10,g11,g12;
|
---|
4146 | g0 = pow(rho,4.0/3.0);
|
---|
4147 | g1 = 1/g0;
|
---|
4148 | g2 = sqrt(gamma);
|
---|
4149 | g3 = asinh(g1*g2);
|
---|
4150 | g4 = 1/a_x;
|
---|
4151 | g5 = 1/pow(rho,4.0/3.0*d);
|
---|
4152 | g6 = d/2.0;
|
---|
4153 | g7 = pow(gamma,g6);
|
---|
4154 | g8 = -(g4*g5*x_d_coef*g7)+6*b*g1*g3*g2+1;
|
---|
4155 | g9 = 1/pow(rho,8.0/3.0);
|
---|
4156 | g10 = pow(gamma,g6-1);
|
---|
4157 | g11 = -beta+b;
|
---|
4158 | g12 = exp(-c*g9*gamma);
|
---|
4159 | dmpw_dg = g0*(((-(1.0/2.0*g4*d*g5*x_d_coef*g10)+(3*b*g9)/sqrt(g9*
|
---|
4160 | gamma+1)+(3*b*g1*g3)/g2)*(-(g11*g9*gamma*g12)-(g5*x_d_coef*g7)+b*g9
|
---|
4161 | *gamma))/pow(g8,2.0)-(((g11*c*gamma*g12)/pow(rho,16.0/3.0)-(g11*g9*
|
---|
4162 | g12)-(1.0/2.0*d*g5*x_d_coef*g10)+b*g9)/g8));
|
---|
4163 | }
|
---|
4164 | }
|
---|
4165 |
|
---|
4166 | void
|
---|
4167 | mPW91XFunctional::point(const PointInputData &id,
|
---|
4168 | PointOutputData &od)
|
---|
4169 | {
|
---|
4170 | od.zero();
|
---|
4171 |
|
---|
4172 | double mpw, dmpw_dr, dmpw_dg;
|
---|
4173 |
|
---|
4174 | // alpha
|
---|
4175 | spin_contrib(id.a, mpw, dmpw_dr, dmpw_dg);
|
---|
4176 | od.energy = mpw;
|
---|
4177 | od.df_drho_a = dmpw_dr;
|
---|
4178 | od.df_dgamma_aa = dmpw_dg;
|
---|
4179 |
|
---|
4180 | // beta
|
---|
4181 | if (spin_polarized_) {
|
---|
4182 | spin_contrib(id.b, mpw, dmpw_dr, dmpw_dg);
|
---|
4183 | od.energy += mpw;
|
---|
4184 | od.df_drho_b = dmpw_dr;
|
---|
4185 | od.df_dgamma_bb = dmpw_dg;
|
---|
4186 | }
|
---|
4187 | else {
|
---|
4188 | od.energy += mpw;
|
---|
4189 | od.df_drho_b = od.df_drho_a;
|
---|
4190 | od.df_dgamma_bb = od.df_dgamma_aa;
|
---|
4191 | }
|
---|
4192 | }
|
---|
4193 |
|
---|
4194 | /////////////////////////////////////////////////////////////////////////////
|
---|
4195 | // Perdew-Wang (PW86) Exchange Functional
|
---|
4196 | // J. P. Perdew and Y. Wang, Phys. Rev. B, 33, 8800, 1986.
|
---|
4197 | //
|
---|
4198 | // Coded by Matt Leininger
|
---|
4199 |
|
---|
4200 | static ClassDesc PW86XFunctional_cd(
|
---|
4201 | typeid(PW86XFunctional),"PW86XFunctional",1,"public DenFunctional",
|
---|
4202 | 0, create<PW86XFunctional>, create<PW86XFunctional>);
|
---|
4203 |
|
---|
4204 | PW86XFunctional::PW86XFunctional(StateIn& s):
|
---|
4205 | SavableState(s),
|
---|
4206 | DenFunctional(s)
|
---|
4207 | {
|
---|
4208 | s.get(m_);
|
---|
4209 | s.get(a_);
|
---|
4210 | s.get(b_);
|
---|
4211 | s.get(c_);
|
---|
4212 | }
|
---|
4213 |
|
---|
4214 | PW86XFunctional::PW86XFunctional()
|
---|
4215 | {
|
---|
4216 | init_constants();
|
---|
4217 | }
|
---|
4218 |
|
---|
4219 | PW86XFunctional::PW86XFunctional(const Ref<KeyVal>& keyval):
|
---|
4220 | DenFunctional(keyval)
|
---|
4221 | {
|
---|
4222 | init_constants();
|
---|
4223 | m_ = keyval->doublevalue("m", KeyValValuedouble(m_));
|
---|
4224 | a_ = keyval->doublevalue("a", KeyValValuedouble(a_));
|
---|
4225 | b_ = keyval->doublevalue("b", KeyValValuedouble(b_));
|
---|
4226 | c_ = keyval->doublevalue("c", KeyValValuedouble(c_));
|
---|
4227 | }
|
---|
4228 |
|
---|
4229 | PW86XFunctional::~PW86XFunctional()
|
---|
4230 | {
|
---|
4231 | }
|
---|
4232 |
|
---|
4233 | void
|
---|
4234 | PW86XFunctional::save_data_state(StateOut& s)
|
---|
4235 | {
|
---|
4236 | DenFunctional::save_data_state(s);
|
---|
4237 | s.put(m_);
|
---|
4238 | s.put(a_);
|
---|
4239 | s.put(b_);
|
---|
4240 | s.put(c_);
|
---|
4241 | }
|
---|
4242 |
|
---|
4243 | void
|
---|
4244 | PW86XFunctional::init_constants()
|
---|
4245 | {
|
---|
4246 | m_ = 1./15.;
|
---|
4247 | a_ = 0.0864/m_;
|
---|
4248 | b_ = 14.;
|
---|
4249 | c_ = 0.2;
|
---|
4250 | }
|
---|
4251 |
|
---|
4252 | int
|
---|
4253 | PW86XFunctional::need_density_gradient()
|
---|
4254 | {
|
---|
4255 | return 1;
|
---|
4256 | }
|
---|
4257 |
|
---|
4258 | void
|
---|
4259 | PW86XFunctional::point(const PointInputData &id,
|
---|
4260 | PointOutputData &od)
|
---|
4261 | {
|
---|
4262 | od.zero();
|
---|
4263 |
|
---|
4264 | double rhoa = 2. * id.a.rho;
|
---|
4265 | double k_fa = pow( (3.*M_PI*M_PI*rhoa), (1./3.) );
|
---|
4266 | double rhoa43 = pow(rhoa, (4./3.));
|
---|
4267 | double rhoa13 = pow(rhoa, (1./3.));
|
---|
4268 | double Ax = -3./4.*pow( (3./M_PI), (1./3.) );
|
---|
4269 | double gamma_aa = 2. * sqrt(id.a.gamma);
|
---|
4270 | double sa;
|
---|
4271 | if (rhoa < MIN_DENSITY) sa = 0.;
|
---|
4272 | else sa = gamma_aa/(2. * k_fa * rhoa);
|
---|
4273 | double sa2 = sa*sa;
|
---|
4274 | double sa3 = sa2*sa;
|
---|
4275 | double sa4 = sa2*sa2;
|
---|
4276 | double sa5 = sa4*sa;
|
---|
4277 | double sa6 = sa5*sa;
|
---|
4278 | double F1a = 1. + a_*sa2 + b_*sa4 + c_*sa6;
|
---|
4279 | double Fxa = pow(F1a, m_);
|
---|
4280 | double fpw86xa = Ax * rhoa43 * Fxa;
|
---|
4281 | double ex = 0.5 * fpw86xa;
|
---|
4282 |
|
---|
4283 | if (compute_potential_) {
|
---|
4284 | double dsa_drhoa;
|
---|
4285 | if (rhoa < MIN_DENSITY) dsa_drhoa = 0.;
|
---|
4286 | else dsa_drhoa = -4.*2.*sa/(3.*rhoa);
|
---|
4287 | double dFxa_drhoa = m_ * pow(F1a, (m_-1.)) * ( 2.*a_*sa + 4.*b_*sa3 + 6.*c_*sa5) * dsa_drhoa;
|
---|
4288 | double dfpw86xa_drhoa = Ax * ( 2.*4./3.*rhoa13*Fxa + rhoa43 * dFxa_drhoa );
|
---|
4289 | od.df_drho_a = 0.5 * dfpw86xa_drhoa;
|
---|
4290 |
|
---|
4291 | double sa_dsa_dgamma_aa;
|
---|
4292 | if (rhoa < MIN_DENSITY) sa_dsa_dgamma_aa = 0.;
|
---|
4293 | else {
|
---|
4294 | sa_dsa_dgamma_aa = 2./pow( (2.*k_fa*rhoa), 2.);
|
---|
4295 | }
|
---|
4296 | double dFxa_dgamma_aa = m_ * pow(F1a, (m_-1.)) *
|
---|
4297 | ( 2.*a_ + 4.*b_*sa2 + 6.*c_*sa4 ) * sa_dsa_dgamma_aa ;
|
---|
4298 | double dfpw86xa_dgamma_aa = Ax * rhoa43 * dFxa_dgamma_aa;
|
---|
4299 | od.df_dgamma_aa = 0.5 * dfpw86xa_dgamma_aa;
|
---|
4300 |
|
---|
4301 | od.df_drho_b = od.df_drho_a;
|
---|
4302 | od.df_dgamma_bb = od.df_dgamma_aa;
|
---|
4303 | od.df_dgamma_ab = 0.;
|
---|
4304 | }
|
---|
4305 |
|
---|
4306 | if (spin_polarized_) {
|
---|
4307 | double rhob = 2. * id.b.rho;
|
---|
4308 | double k_fb = pow( (3.*M_PI*M_PI*rhob), (1./3.) );
|
---|
4309 | double rhob43 = pow(rhob, (4./3.));
|
---|
4310 | double rhob13 = pow(rhob, (1./3.));
|
---|
4311 | double gamma_bb = 2.*sqrt(id.b.gamma);
|
---|
4312 | double sb;
|
---|
4313 | if (rhob < MIN_DENSITY) sb = 0.;
|
---|
4314 | else sb = gamma_bb/(2.*k_fb*rhob);
|
---|
4315 | double sb2 = sb*sb;
|
---|
4316 | double sb3 = sb2*sb;
|
---|
4317 | double sb4 = sb3*sb;
|
---|
4318 | double sb5 = sb4*sb;
|
---|
4319 | double sb6 = sb5*sb;
|
---|
4320 | double F1b = 1. + a_*sb2 + b_*sb4 + c_*sb6;
|
---|
4321 | double Fxb = pow(F1b, m_);
|
---|
4322 | double fpw86xb = Ax * rhob43 * Fxb;
|
---|
4323 | ex += 0.5 * fpw86xb;
|
---|
4324 |
|
---|
4325 | if (compute_potential_) {
|
---|
4326 | double dsb_drhob;
|
---|
4327 | if (rhob < MIN_DENSITY) dsb_drhob = 0.;
|
---|
4328 | else dsb_drhob = -4.*2.*sb/(3.*rhob);
|
---|
4329 | double dFxb_drhob = m_ * pow(F1b, (m_-1.))
|
---|
4330 | * (2.*a_*sb + 4.*b_*sb3 + 6.*c_*sb5) * dsb_drhob;
|
---|
4331 | double dfpw86xb = Ax * ( 2.*4./3.*rhob13*Fxb + rhob43 * dFxb_drhob );
|
---|
4332 | od.df_drho_b = 0.5 * dfpw86xb;
|
---|
4333 |
|
---|
4334 | double sb_dsb_dgamma_bb;
|
---|
4335 | if (rhob < MIN_DENSITY) sb_dsb_dgamma_bb = 0.;
|
---|
4336 | else {
|
---|
4337 | sb_dsb_dgamma_bb = 2./pow( (2.*k_fb*rhob), 2.);
|
---|
4338 | }
|
---|
4339 | double dFxb_dgamma_bb = m_ * pow(F1b, (m_-1.)) *
|
---|
4340 | ( 2.*a_ + 4.*b_*sb2 + 6.*c_*sb4) * sb_dsb_dgamma_bb;
|
---|
4341 | double dfpw86xb_dgamma_bb = Ax * rhob43 * dFxb_dgamma_bb;
|
---|
4342 | od.df_dgamma_bb = 0.5 * dfpw86xb_dgamma_bb;
|
---|
4343 | }
|
---|
4344 | }
|
---|
4345 | else ex += ex;
|
---|
4346 |
|
---|
4347 | od.energy = ex;
|
---|
4348 | }
|
---|
4349 |
|
---|
4350 | /////////////////////////////////////////////////////////////////////////////
|
---|
4351 | // Gill 1996 (G96) Exchange Functional
|
---|
4352 | // P. M. W. Gill, Mol. Phys. 89, 433, 1996.
|
---|
4353 | //
|
---|
4354 | // Coded by Matt Leininger
|
---|
4355 |
|
---|
4356 | static ClassDesc G96XFunctional_cd(
|
---|
4357 | typeid(G96XFunctional),"G96XFunctional",1,"public DenFunctional",
|
---|
4358 | 0, create<G96XFunctional>, create<G96XFunctional>);
|
---|
4359 |
|
---|
4360 | G96XFunctional::G96XFunctional(StateIn& s):
|
---|
4361 | SavableState(s),
|
---|
4362 | DenFunctional(s)
|
---|
4363 | {
|
---|
4364 | s.get(b_);
|
---|
4365 | }
|
---|
4366 |
|
---|
4367 | G96XFunctional::G96XFunctional()
|
---|
4368 | {
|
---|
4369 | init_constants();
|
---|
4370 | }
|
---|
4371 |
|
---|
4372 | G96XFunctional::G96XFunctional(const Ref<KeyVal>& keyval):
|
---|
4373 | DenFunctional(keyval)
|
---|
4374 | {
|
---|
4375 | init_constants();
|
---|
4376 | b_ = keyval->doublevalue("b", KeyValValuedouble(b_));
|
---|
4377 | }
|
---|
4378 |
|
---|
4379 | G96XFunctional::~G96XFunctional()
|
---|
4380 | {
|
---|
4381 | }
|
---|
4382 |
|
---|
4383 | void
|
---|
4384 | G96XFunctional::save_data_state(StateOut& s)
|
---|
4385 | {
|
---|
4386 | DenFunctional::save_data_state(s);
|
---|
4387 | s.put(b_);
|
---|
4388 | }
|
---|
4389 |
|
---|
4390 | void
|
---|
4391 | G96XFunctional::init_constants()
|
---|
4392 | {
|
---|
4393 | b_ = 1./137.;
|
---|
4394 | }
|
---|
4395 |
|
---|
4396 | int
|
---|
4397 | G96XFunctional::need_density_gradient()
|
---|
4398 | {
|
---|
4399 | return 1;
|
---|
4400 | }
|
---|
4401 |
|
---|
4402 | void
|
---|
4403 | G96XFunctional::point(const PointInputData &id,
|
---|
4404 | PointOutputData &od)
|
---|
4405 | {
|
---|
4406 | od.zero();
|
---|
4407 |
|
---|
4408 | double rhoa = id.a.rho;
|
---|
4409 | double rhoa43 = pow(rhoa, (4./3.));
|
---|
4410 | double rhoa13 = pow(rhoa, (1./3.));
|
---|
4411 | double gamma_aa = sqrt(id.a.gamma);
|
---|
4412 | double gamma_aa32 = pow(gamma_aa, 1.5);
|
---|
4413 | double alpha = -1.5 * pow( (3./(4.*M_PI)), (1./3.) );
|
---|
4414 | double gg96a;
|
---|
4415 | double fxg96a;
|
---|
4416 | if (rhoa < MIN_DENSITY) gg96a = fxg96a = 0.;
|
---|
4417 | else {
|
---|
4418 | gg96a = alpha - b_*gamma_aa32/(rhoa*rhoa);
|
---|
4419 | fxg96a = rhoa43 * gg96a;
|
---|
4420 | }
|
---|
4421 | double ex = fxg96a;
|
---|
4422 |
|
---|
4423 | if (compute_potential_) {
|
---|
4424 |
|
---|
4425 | double dfxg96a_drhoa;
|
---|
4426 | if (rhoa < MIN_DENSITY) dfxg96a_drhoa = 0.;
|
---|
4427 | else {
|
---|
4428 | double dgg96a_drhoa = 2.*b_*gamma_aa32/(rhoa*rhoa*rhoa);
|
---|
4429 | dfxg96a_drhoa = 4./3.*rhoa13*gg96a + rhoa43*dgg96a_drhoa;
|
---|
4430 | }
|
---|
4431 | od.df_drho_a = dfxg96a_drhoa;
|
---|
4432 | od.df_drho_b = od.df_drho_a;
|
---|
4433 |
|
---|
4434 | double dfxg96a_dgamma_aa;
|
---|
4435 | // The derivative of the G96X functional with respect to gamma_aa or bb
|
---|
4436 | // as implemented should go to infinity as gamma goes to zero.
|
---|
4437 | // However, the derivative gamma terms are eventually contracted with quantities
|
---|
4438 | // that have a sqrt(gamma) in the numerator and therefore the overall limit
|
---|
4439 | // is zero.
|
---|
4440 | if (gamma_aa < MIN_GAMMA || rhoa < MIN_DENSITY ) dfxg96a_dgamma_aa = 0.;
|
---|
4441 | else {
|
---|
4442 | double dgg96a_dgamma_aa = -3.*b_ / ( 4.*rhoa*rhoa*sqrt(gamma_aa) );
|
---|
4443 | dfxg96a_dgamma_aa = rhoa43 * dgg96a_dgamma_aa;
|
---|
4444 | }
|
---|
4445 | od.df_dgamma_aa = dfxg96a_dgamma_aa;
|
---|
4446 | od.df_dgamma_bb = od.df_dgamma_aa;
|
---|
4447 | od.df_dgamma_ab = 0.;
|
---|
4448 | }
|
---|
4449 |
|
---|
4450 | if (spin_polarized_) {
|
---|
4451 | double rhob = id.b.rho;
|
---|
4452 | double rhob43 = pow(rhob, (4./3.));
|
---|
4453 | double rhob13 = pow(rhob, (1./3.));
|
---|
4454 | double gamma_bb = sqrt(id.b.gamma);
|
---|
4455 | double gamma_bb32 = pow(gamma_bb, 1.5);
|
---|
4456 | double gg96b;
|
---|
4457 | double fxg96b;
|
---|
4458 | if (rhob < MIN_DENSITY) gg96b = fxg96b = 0.;
|
---|
4459 | else {
|
---|
4460 | gg96b = alpha - b_*gamma_bb32/(rhob*rhob);
|
---|
4461 | fxg96b = rhob43 * gg96b;
|
---|
4462 | }
|
---|
4463 | ex += fxg96b;
|
---|
4464 |
|
---|
4465 | if (compute_potential_) {
|
---|
4466 | double dfxg96b_drhob;
|
---|
4467 | if (rhob < MIN_DENSITY) dfxg96b_drhob = 0.;
|
---|
4468 | else {
|
---|
4469 | double dgg96b_drhob = 2.*b_*gamma_bb32/(rhob*rhob*rhob);
|
---|
4470 | dfxg96b_drhob = 4./3.*rhob13*gg96b + rhob43*dgg96b_drhob;
|
---|
4471 | }
|
---|
4472 | od.df_drho_b = dfxg96b_drhob;
|
---|
4473 |
|
---|
4474 | double dfxg96b_dgamma_bb;
|
---|
4475 | // See comment above with regard to correct limits.
|
---|
4476 | if (gamma_bb < MIN_GAMMA || rhob < MIN_DENSITY) dfxg96b_dgamma_bb=0.;
|
---|
4477 | else {
|
---|
4478 | double dgg96b_dgamma_bb = -3.*b_ / (4.*rhob*rhob*sqrt(gamma_bb));
|
---|
4479 | dfxg96b_dgamma_bb = rhob43 * dgg96b_dgamma_bb;
|
---|
4480 | }
|
---|
4481 | od.df_dgamma_bb = dfxg96b_dgamma_bb;
|
---|
4482 | }
|
---|
4483 | }
|
---|
4484 | else ex += ex;
|
---|
4485 |
|
---|
4486 |
|
---|
4487 | od.energy = ex;
|
---|
4488 | }
|
---|
4489 |
|
---|
4490 | /////////////////////////////////////////////////////////////////////////////
|
---|
4491 |
|
---|
4492 | // Local Variables:
|
---|
4493 | // mode: c++
|
---|
4494 | // c-file-style: "CLJ"
|
---|
4495 | // End:
|
---|