source: src/Jobs/mpqc_extract.cc@ 9a4772

Action_Thermostats Add_SelectAtomByNameAction Adding_Graph_to_ChangeBondActions Adding_MD_integration_tests Adding_StructOpt_integration_tests AutomationFragmentation_failures Candidate_v1.6.1 Candidate_v1.7.0 ChemicalSpaceEvaluator EmpiricalPotential_contain_HomologyGraph_documentation Enhance_userguide Enhanced_StructuralOptimization Enhanced_StructuralOptimization_continued Example_ManyWaysToTranslateAtom Exclude_Hydrogens_annealWithBondGraph Fix_Verbose_Codepatterns ForceAnnealing_oldresults ForceAnnealing_with_BondGraph ForceAnnealing_with_BondGraph_continued ForceAnnealing_with_BondGraph_continued_betteresults ForceAnnealing_with_BondGraph_contraction-expansion Gui_displays_atomic_force_velocity IndependentFragmentGrids_IntegrationTest JobMarket_RobustOnKillsSegFaults JobMarket_StableWorkerPool PythonUI_with_named_parameters Recreated_GuiChecks StoppableMakroAction TremoloParser_IncreasedPrecision TremoloParser_MultipleTimesteps stable
Last change on this file since 9a4772 was fbf005, checked in by Frederik Heber <heber@…>, 9 years ago

HUGE: Extracted libmolecuilder_mpqc that is now linked to poolworker.

  • This stops the problems with MemDebug and mpqc when linking against libLinearAlgebra in debug mode: static global variables in libLinAlg are allocated using MemDebug (prefixed with a checksum) but are deallocated using normal libc's free() on exit. This causes an invalid free() as the ptr given to free points inside the block and not at its start. The problem comes from mpqc's use of own new and delete implementation. I'm guessing its reference counting is the culprit. Hence, we need to separate these two compilations from another in different units/libraries. Therefore, we have split off libmolecuilder_mpqc, .._mpqc_extract and additionally place the MPQCJob::Work() implementation inside libMolecuilderJobs_Work.
  • libmolecuilder_mpqc contains all MPQC's code in mpqc.cc (and linked libraries) that is not the main() function.
  • libmolecuilder_mpqc_extract contains functions that extract results such as energies, forces, charge grids from the obtained mpqc solution. These were added by myself to the mpqc code before.
  • molecuilder_mpqc is then linked against a NoOp .._extract library version. Thereby, it does not use any of the Molecuilder or related libraries and does not come in contact with MemDebug.
  • molecuilder_poolworker however is linked with the full .._extract library and hence performs these extractions.
  • poolworker now executes MPQCJob, MPQCCommandJob, and VMGJob and therefore needs to enforce binding to all of them.
  • TESTS: renamed molecuilder_mpqc.in to molecuilder_poolworker.in.
  • Property mode set to 100644
File size: 21.9 KB
Line 
1//
2// mpqc_extract.cc
3//
4// Copyright (C) 1996 Limit Point Systems, Inc.
5//
6// Author: Edward Seidl <seidl@janed.com>
7// Maintainer: LPS
8//
9// This file is part of MPQC.
10//
11// MPQC is free software; you can redistribute it and/or modify
12// it under the terms of the GNU General Public License as published by
13// the Free Software Foundation; either version 2, or (at your option)
14// any later version.
15//
16// MPQC 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 General Public License for more details.
20//
21// You should have received a copy of the GNU General Public License
22// along with the MPQC; see the file COPYING. 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// \note This was extracted from \file mpqc.cc for refactoring into a library.
28
29#ifdef HAVE_CONFIG_H
30#include <scconfig.h>
31#endif
32
33#ifdef HAVE_JOBMARKET
34// include headers that implement a archive in simple text format
35// otherwise BOOST_CLASS_EXPORT_IMPLEMENT has no effect
36#include <boost/archive/text_oarchive.hpp>
37#include <boost/archive/text_iarchive.hpp>
38
39#include "JobMarket/Results/FragmentResult.hpp"
40#include "JobMarket/poolworker_main.hpp"
41
42#include "chemistry/qc/scf/scfops.h"
43
44#ifdef HAVE_MPQCDATA
45#include "Jobs/MPQCJob.hpp"
46#include "Fragmentation/Summation/Containers/MPQCData.hpp"
47
48#include <chemistry/qc/basis/obint.h>
49#include <chemistry/qc/basis/symmint.h>
50#endif
51
52#include <algorithm>
53#include <stdlib.h>
54#endif
55
56#include <chemistry/qc/scf/linkage.h>
57#include <chemistry/qc/dft/linkage.h>
58#include <chemistry/qc/mbpt/linkage.h>
59#ifdef HAVE_SC_SRC_LIB_CHEMISTRY_QC_MBPTR12
60# include <chemistry/qc/mbptr12/linkage.h>
61#endif
62#ifdef HAVE_SC_SRC_LIB_CHEMISTRY_QC_CINTS
63# include <chemistry/qc/cints/linkage.h>
64#endif
65//#include <chemistry/qc/psi/linkage.h>
66#include <util/state/linkage.h>
67#ifdef HAVE_SC_SRC_LIB_CHEMISTRY_QC_CC
68# include <chemistry/qc/cc/linkage.h>
69#endif
70#ifdef HAVE_SC_SRC_LIB_CHEMISTRY_QC_PSI
71# include <chemistry/qc/psi/linkage.h>
72#endif
73#ifdef HAVE_SC_SRC_LIB_CHEMISTRY_QC_INTCCA
74# include <chemistry/qc/intcca/linkage.h>
75#endif
76
77#include "mpqc_extract.h"
78
79using namespace std;
80using namespace sc;
81
82static int getCoreElectrons(const int z)
83{
84 int n=0;
85 if (z > 2) n += 2;
86 if (z > 10) n += 8;
87 if (z > 18) n += 8;
88 if (z > 30) n += 10;
89 if (z > 36) n += 8;
90 if (z > 48) n += 10;
91 if (z > 54) n += 8;
92 return n;
93}
94
95/** Finds the region index to a given timer region name.
96 *
97 * @param nregion number of regions
98 * @param region_names array with name of each region
99 * @param name name of desired region
100 * @return index of desired region in array
101 */
102int findTimerRegion(const int &nregion, const char **&region_names, const char *name)
103{
104 int region=0;
105 for (;region<nregion;++region) {
106 //std::cout << "Comparing " << region_names[region] << " and " << name << "." << std::endl;
107 if (strcmp(region_names[region], name) == 0)
108 break;
109 }
110 if (region == nregion)
111 region = 0;
112 return region;
113}
114
115/** Extractor function that is called after all calculations have been made.
116 *
117 * \param data result structure to fill
118 */
119void extractResults(
120 Ref<MolecularEnergy> &mole,
121 void *_data
122 )
123{
124 MPQCData &data = *static_cast<MPQCData *>(_data);
125 Ref<Wavefunction> wfn;
126 wfn << mole;
127// ExEnv::out0() << "The number of atomic orbitals: " << wfn->ao_dimension()->n() << endl;
128// ExEnv::out0() << "The AO density matrix is ";
129// wfn->ao_density()->print(ExEnv::out0());
130// ExEnv::out0() << "The natural density matrix is ";
131// wfn->natural_density()->print(ExEnv::out0());
132// ExEnv::out0() << "The Gaussian basis is " << wfn->basis()->name() << endl;
133// ExEnv::out0() << "The Gaussians sit at the following centers: " << endl;
134// for (int nr = 0; nr< wfn->basis()->ncenter(); ++nr) {
135// ExEnv::out0() << nr << " basis function has its center at ";
136// for (int i=0; i < 3; ++i)
137// ExEnv::out0() << wfn->basis()->r(nr,i) << "\t";
138// ExEnv::out0() << endl;
139// }
140 // store accuracies
141 data.accuracy = mole->value_result().actual_accuracy();
142 data.desired_accuracy = mole->value_result().desired_accuracy();
143 // print the energy
144 data.energies.total = wfn->energy();
145 data.energies.nuclear_repulsion = wfn->nuclear_repulsion_energy();
146 {
147 CLHF *clhf = dynamic_cast<CLHF*>(wfn.pointer());
148 if (clhf != NULL) {
149 double ex, ec;
150 clhf->two_body_energy(ec, ex);
151 data.energies.electron_coulomb = ec;
152 data.energies.electron_exchange = ex;
153 clhf = NULL;
154 } else {
155 ExEnv::out0() << "INFO: There is no direct CLHF information available." << endl;
156 data.energies.electron_coulomb = 0.;
157 data.energies.electron_exchange = 0.;
158 }
159 }
160 SCF *scf = NULL;
161 {
162 MBPT2 *mbpt2 = dynamic_cast<MBPT2*>(wfn.pointer());
163 if (mbpt2 != NULL) {
164 data.energies.correlation = mbpt2->corr_energy();
165 scf = mbpt2->ref().pointer();
166 CLHF *clhf = dynamic_cast<CLHF*>(scf);
167 if (clhf != NULL) {
168 double ex, ec;
169 clhf->two_body_energy(ec, ex);
170 data.energies.electron_coulomb = ec;
171 data.energies.electron_exchange = ex;
172 clhf = NULL;
173 } else {
174 ExEnv::out0() << "INFO: There is no reference CLHF information available either." << endl;
175 data.energies.electron_coulomb = 0.;
176 data.energies.electron_exchange = 0.;
177 }
178 mbpt2 = 0;
179 } else {
180 ExEnv::out0() << "INFO: There is no MBPT2 information available." << endl;
181 data.energies.correlation = 0.;
182 scf = dynamic_cast<SCF*>(wfn.pointer());
183 if (scf == NULL)
184 abort();
185 }
186 }
187 {
188 // taken from clscf.cc: CLSCF::scf_energy() (but see also Szabo/Ostlund)
189
190 RefSymmSCMatrix t = scf->overlap();
191 RefSymmSCMatrix cl_dens_ = scf->ao_density();
192
193 SCFEnergy *eop = new SCFEnergy;
194 eop->reference();
195 if (t.dim()->equiv(cl_dens_.dim())) {
196 Ref<SCElementOp2> op = eop;
197 t.element_op(op,cl_dens_);
198 op=0;
199 }
200 eop->dereference();
201
202 data.energies.overlap = eop->result();
203
204 delete eop;
205 t = 0;
206 cl_dens_ = 0;
207 }
208 {
209 // taken from Wavefunction::core_hamiltonian()
210 RefSymmSCMatrix hao(scf->basis()->basisdim(), scf->basis()->matrixkit());
211 hao.assign(0.0);
212 Ref<PetiteList> pl = scf->integral()->petite_list();
213 Ref<SCElementOp> hc =
214 new OneBodyIntOp(new SymmOneBodyIntIter(scf->integral()->kinetic(), pl));
215 hao.element_op(hc);
216 hc=0;
217
218 RefSymmSCMatrix h(scf->so_dimension(), scf->basis_matrixkit());
219 pl->symmetrize(hao,h);
220
221 // taken from clscf.cc: CLSCF::scf_energy() (but see also Szabo/Ostlund)
222 RefSymmSCMatrix cl_dens_ = scf->ao_density();
223
224 SCFEnergy *eop = new SCFEnergy;
225 eop->reference();
226 if (h.dim()->equiv(cl_dens_.dim())) {
227 Ref<SCElementOp2> op = eop;
228 h.element_op(op,cl_dens_);
229 op=0;
230 }
231 eop->dereference();
232
233 data.energies.kinetic = 2.*eop->result();
234
235 delete eop;
236 hao = 0;
237 h = 0;
238 cl_dens_ = 0;
239 }
240 {
241 // set to potential energy between nuclei and electron charge distribution
242 RefSymmSCMatrix hao(scf->basis()->basisdim(), scf->basis()->matrixkit());
243 hao.assign(0.0);
244 Ref<PetiteList> pl = scf->integral()->petite_list();
245 Ref<SCElementOp> hc =
246 new OneBodyIntOp(new SymmOneBodyIntIter(scf->integral()->nuclear(), pl));
247 hao.element_op(hc);
248 hc=0;
249
250 RefSymmSCMatrix h(scf->so_dimension(), scf->basis_matrixkit());
251 pl->symmetrize(hao,h);
252
253 // taken from clscf.cc: CLSCF::scf_energy() (but see also Szabo/Ostlund)
254 RefSymmSCMatrix cl_dens_ = scf->ao_density();
255
256 SCFEnergy *eop = new SCFEnergy;
257 eop->reference();
258 if (h.dim()->equiv(cl_dens_.dim())) {
259 Ref<SCElementOp2> op = eop;
260 h.element_op(op,cl_dens_);
261 op=0;
262 }
263 eop->dereference();
264
265 data.energies.hcore = 2.*eop->result();
266
267 delete eop;
268 hao = 0;
269 h = 0;
270 cl_dens_ = 0;
271 }
272 ExEnv::out0() << "total is " << data.energies.total << endl;
273 ExEnv::out0() << "nuclear_repulsion is " << data.energies.nuclear_repulsion << endl;
274 ExEnv::out0() << "electron_coulomb is " << data.energies.electron_coulomb << endl;
275 ExEnv::out0() << "electron_exchange is " << data.energies.electron_exchange << endl;
276 ExEnv::out0() << "correlation is " << data.energies.correlation << endl;
277 ExEnv::out0() << "overlap is " << data.energies.overlap << endl;
278 ExEnv::out0() << "kinetic is " << data.energies.kinetic << endl;
279 ExEnv::out0() << "hcore is " << data.energies.hcore << endl;
280 ExEnv::out0() << "sum is " <<
281 data.energies.nuclear_repulsion
282 + data.energies.electron_coulomb
283 + data.energies.electron_exchange
284 + data.energies.correlation
285 + data.energies.kinetic
286 + data.energies.hcore
287 << endl;
288
289 ExEnv::out0() << endl << indent
290 << scprintf("Value of the MolecularEnergy: %15.10f",
291 mole->energy())
292 << endl;
293 // print the gradient
294 RefSCVector grad;
295 if (mole->gradient_result().computed()) {
296 grad = mole->gradient_result().result_noupdate();
297 }
298 // gradient calculation needs to be activated in the configuration
299 // some methods such as open shell MBPT2 do not allow for gradient calc.
300// else {
301// grad = mole->gradient();
302// }
303 if (grad.nonnull()) {
304 data.forces.resize(grad.dim()/3);
305 for (int j=0;j<grad.dim()/3; ++j) {
306 data.forces[j].resize(3, 0.);
307 }
308 ExEnv::out0() << "Gradient of the MolecularEnergy:" << std::endl;
309 for (int j=0;j<grad.dim()/3; ++j) {
310 ExEnv::out0() << "\t";
311 for (int i=0; i< 3; ++i) {
312 data.forces[j][i] = grad[3*j+i];
313 ExEnv::out0() << grad[3*j+i] << "\t";
314 }
315 ExEnv::out0() << endl;
316 }
317 } else {
318 ExEnv::out0() << "INFO: There is no gradient information available." << endl;
319 }
320 grad = NULL;
321
322 {
323 // eigenvalues (this only works if we have a OneBodyWavefunction, i.e. SCF procedure)
324 // eigenvalues seem to be invalid for unrestricted SCF calculation
325 // (see UnrestrictedSCF::eigenvalues() implementation)
326 UnrestrictedSCF *uscf = dynamic_cast<UnrestrictedSCF*>(wfn.pointer());
327 if ((scf != NULL) && (uscf == NULL)) {
328// const double scfernergy = scf->energy();
329 RefDiagSCMatrix evals = scf->eigenvalues();
330
331 ExEnv::out0() << "Eigenvalues:" << endl;
332 for(int i=0;i<wfn->oso_dimension(); ++i) {
333 data.energies.eigenvalues.push_back(evals(i));
334 ExEnv::out0() << i << "th eigenvalue is " << evals(i) << endl;
335 }
336 } else {
337 ExEnv::out0() << "INFO: There is no eigenvalue information available." << endl;
338 }
339 }
340 // we do sample the density only on request
341 {
342 // fill positions and charges (NO LONGER converting from bohr radii to angstroem)
343 const double AtomicLengthToAngstroem = 1.;//0.52917721;
344 data.positions.reserve(wfn->molecule()->natom());
345 data.atomicnumbers.reserve(wfn->molecule()->natom());
346 data.charges.reserve(wfn->molecule()->natom());
347 for (int iatom=0;iatom < wfn->molecule()->natom(); ++iatom) {
348 data.atomicnumbers.push_back(wfn->molecule()->Z(iatom));
349 double charge = wfn->molecule()->Z(iatom);
350 if (data.DoValenceOnly == MPQCData::DoSampleValenceOnly)
351 charge -= getCoreElectrons((int)charge);
352 data.charges.push_back(charge);
353 std::vector<double> pos(3, 0.);
354 for (int j=0;j<3;++j)
355 pos[j] = wfn->molecule()->r(iatom, j)*AtomicLengthToAngstroem;
356 data.positions.push_back(pos);
357 }
358 ExEnv::out0() << "We have "
359 << data.positions.size() << " positions and "
360 << data.charges.size() << " charges." << endl;
361 }
362 if (data.DoLongrange) {
363 if (data.sampled_grid.level != 0)
364 {
365 // we now need to sample the density on the grid
366 // 1. get max and min over all basis function positions
367 assert( scf->basis()->ncenter() > 0 );
368 SCVector3 bmin( scf->basis()->r(0,0), scf->basis()->r(0,1), scf->basis()->r(0,2) );
369 SCVector3 bmax( scf->basis()->r(0,0), scf->basis()->r(0,1), scf->basis()->r(0,2) );
370 for (int nr = 1; nr< scf->basis()->ncenter(); ++nr) {
371 for (int i=0; i < 3; ++i) {
372 if (scf->basis()->r(nr,i) < bmin(i))
373 bmin(i) = scf->basis()->r(nr,i);
374 if (scf->basis()->r(nr,i) > bmax(i))
375 bmax(i) = scf->basis()->r(nr,i);
376 }
377 }
378 ExEnv::out0() << "Basis min is at " << bmin << " and max is at " << bmax << endl;
379
380 // 2. choose an appropriately large grid
381 // we have to pay attention to capture the right amount of the exponential decay
382 // and also to have a power of two size of the grid at best
383 SCVector3 boundaryV(5.); // boundary extent around compact domain containing basis functions
384 bmin -= boundaryV;
385 bmax += boundaryV;
386 for (size_t i=0;i<3;++i) {
387 if (bmin(i) < data.sampled_grid.begin[i])
388 bmin(i) = data.sampled_grid.begin[i];
389 if (bmax(i) > data.sampled_grid.end[i])
390 bmax(i) = data.sampled_grid.end[i];
391 }
392 // set the non-zero window of the sampled_grid
393 data.sampled_grid.setWindow(bmin.data(), bmax.data());
394
395 // for the moment we always generate a grid of full size
396 // (NO LONGER converting grid dimensions from angstroem to bohr radii)
397 const double AtomicLengthToAngstroem = 1.;//0.52917721;
398 SCVector3 min;
399 SCVector3 max;
400 SCVector3 delta;
401 size_t samplepoints[3];
402 // due to periodic boundary conditions, we don't need gridpoints-1 here
403 // TODO: in case of open boundary conditions, we need data on the right
404 // hand side boundary as well
405// const int gridpoints = data.sampled_grid.getGridPointsPerAxis();
406 for (size_t i=0;i<3;++i) {
407 min(i) = data.sampled_grid.begin_window[i]/AtomicLengthToAngstroem;
408 max(i) = data.sampled_grid.end_window[i]/AtomicLengthToAngstroem;
409 delta(i) = data.sampled_grid.getDeltaPerAxis(i)/AtomicLengthToAngstroem;
410 samplepoints[i] = data.sampled_grid.getWindowGridPointsPerAxis(i);
411 }
412 ExEnv::out0() << "Grid starts at " << min
413 << " and ends at " << max
414 << " with a delta of " << delta
415 << " to get "
416 << samplepoints[0] << ","
417 << samplepoints[1] << ","
418 << samplepoints[2] << " samplepoints."
419 << endl;
420 assert( data.sampled_grid.sampled_grid.size() == samplepoints[0]*samplepoints[1]*samplepoints[2]);
421
422 // 3. sample the atomic density
423 const double element_volume_conversion =
424 1./AtomicLengthToAngstroem/AtomicLengthToAngstroem/AtomicLengthToAngstroem;
425 SCVector3 r = min;
426
427 std::set<int> valence_indices;
428 RefDiagSCMatrix evals = scf->eigenvalues();
429 if (data.DoValenceOnly == MPQCData::DoSampleValenceOnly) {
430 // find valence orbitals
431// std::cout << "All Eigenvalues:" << std::endl;
432// for(int i=0;i<wfn->oso_dimension(); ++i)
433// std::cout << i << "th eigenvalue is " << evals(i) << std::endl;
434// int n_electrons = scf->nelectron();
435 int n_core_electrons = wfn->molecule()->n_core_electrons();
436 std::set<double> evals_sorted;
437 {
438 int i=0;
439 double first_positive_ev = std::numeric_limits<double>::max();
440 for(i=0;i<wfn->oso_dimension(); ++i) {
441 if (evals(i) < 0.)
442 evals_sorted.insert(evals(i));
443 else
444 first_positive_ev = std::min(first_positive_ev, (double)evals(i));
445 }
446 // add the first positive for the distance
447 evals_sorted.insert(first_positive_ev);
448 }
449 std::set<double> evals_distances;
450 std::set<double>::const_iterator advancer = evals_sorted.begin();
451 std::set<double>::const_iterator iter = advancer++;
452 for(;advancer != evals_sorted.end(); ++advancer,++iter)
453 evals_distances.insert((*advancer)-(*iter));
454 const double largest_distance = *(evals_distances.rbegin());
455 ExEnv::out0() << "Largest distance between EV is " << largest_distance << std::endl;
456 advancer = evals_sorted.begin();
457 iter = advancer++;
458 for(;advancer != evals_sorted.begin(); ++advancer,++iter)
459 if (fabs(fabs((*advancer)-(*iter)) - largest_distance) < 1e-10)
460 break;
461 assert( advancer != evals_sorted.begin() );
462 const double last_core_ev = (*iter);
463 ExEnv::out0() << "Last core EV might be " << last_core_ev << std::endl;
464 ExEnv::out0() << "First valence index is " << n_core_electrons/2 << std::endl;
465 for(int i=n_core_electrons/2;i<wfn->oso_dimension(); ++i)
466 if (evals(i) > last_core_ev)
467 valence_indices.insert(i);
468// {
469// int i=0;
470// std::cout << "Valence eigenvalues:" << std::endl;
471// for (std::set<int>::const_iterator iter = valence_indices.begin();
472// iter != valence_indices.end(); ++iter)
473// std::cout << i++ << "th eigenvalue is " << (*iter) << std::endl;
474// }
475 } else {
476 // just insert all indices
477 for(int i=0;i<wfn->oso_dimension(); ++i)
478 valence_indices.insert(i);
479 }
480
481 // testing alternative routine from SCF::so_density()
482 RefSCMatrix oso_vector = scf->oso_eigenvectors();
483 RefSCMatrix vector = scf->so_to_orthog_so().t() * oso_vector;
484 oso_vector = 0;
485 RefSymmSCMatrix occ(scf->oso_dimension(), scf->basis_matrixkit());
486 occ.assign(0.0);
487 for (std::set<int>::const_iterator iter = valence_indices.begin();
488 iter != valence_indices.end(); ++iter) {
489 const int i = *iter;
490 occ(i,i) = scf->occupation(i);
491 ExEnv::out0() << "# " << i << " has ev of " << evals(i) << ", occupied by " << scf->occupation(i) << std::endl;
492 }
493 RefSymmSCMatrix d2(scf->so_dimension(), scf->basis_matrixkit());
494 d2.assign(0.0);
495 d2.accumulate_transform(vector, occ);
496
497 // taken from scf::density()
498 RefSCMatrix nos
499 = scf->integral()->petite_list()->evecs_to_AO_basis(scf->natural_orbitals());
500 RefDiagSCMatrix nd = scf->natural_density();
501 GaussianBasisSet::ValueData *valdat
502 = new GaussianBasisSet::ValueData(scf->basis(), scf->integral());
503 std::vector<double>::iterator griditer = data.sampled_grid.sampled_grid.begin();
504 const int nbasis = scf->basis()->nbasis();
505 double *bs_values = new double[nbasis];
506
507 // TODO: need to take care when we have periodic boundary conditions.
508 for (size_t x = 0; x < samplepoints[0]; ++x, r.x() += delta(0)) {
509 std::cout << "Sampling now for x=" << r.x() << std::endl;
510 for (size_t y = 0; y < samplepoints[1]; ++y, r.y() += delta(1)) {
511 for (size_t z = 0; z < samplepoints[2]; ++z, r.z() += delta(2)) {
512 scf->basis()->values(r,valdat,bs_values);
513
514 // loop over natural orbitals adding contributions to elec_density
515 double elec_density=0.0;
516 for (int i=0; i<nbasis; ++i) {
517 double tmp = 0.0;
518 for (int j=0; j<nbasis; ++j) {
519 tmp += d2(j,i)*bs_values[j]*bs_values[i];
520 }
521 elec_density += tmp;
522 }
523 const double dens_at_r = elec_density * element_volume_conversion;
524// const double dens_at_r = scf->density(r) * element_volume_conversion;
525
526// if (fabs(dens_at_r) > 1e-4)
527// std::cout << "Electron density at " << r << " is " << dens_at_r << std::endl;
528 if (griditer != data.sampled_grid.sampled_grid.end())
529 *griditer++ = dens_at_r;
530 else
531 std::cerr << "PAST RANGE!" << std::endl;
532 }
533 r.z() = min.z();
534 }
535 r.y() = min.y();
536 }
537 delete[] bs_values;
538 delete valdat;
539 assert( griditer == data.sampled_grid.sampled_grid.end());
540 // normalization of electron charge to equal electron number
541 {
542 double integral_value = 0.;
543 const double volume_element = pow(AtomicLengthToAngstroem,3)*delta(0)*delta(1)*delta(2);
544 for (std::vector<double>::const_iterator diter = data.sampled_grid.sampled_grid.begin();
545 diter != data.sampled_grid.sampled_grid.end(); ++diter)
546 integral_value += *diter;
547 integral_value *= volume_element;
548 int n_electrons = scf->nelectron();
549 if (data.DoValenceOnly == MPQCData::DoSampleValenceOnly)
550 n_electrons -= wfn->molecule()->n_core_electrons();
551 const double normalization =
552 ((integral_value == 0) || (n_electrons == 0)) ?
553 1. : n_electrons/integral_value;
554 std::cout << "Created " << data.sampled_grid.sampled_grid.size() << " grid points"
555 << " with integral value of " << integral_value
556 << " against " << ((data.DoValenceOnly == MPQCData::DoSampleValenceOnly) ? "n_valence_electrons" : "n_electrons")
557 << " of " << n_electrons << "." << std::endl;
558 // with normalization we also get the charge right : -1.
559 for (std::vector<double>::iterator diter = data.sampled_grid.sampled_grid.begin();
560 diter != data.sampled_grid.sampled_grid.end(); ++diter)
561 *diter *= -1.*normalization;
562 }
563 }
564 }
565 scf = 0;
566}
567
568void extractTimings(
569 Ref<RegionTimer> &tim,
570 void *_data)
571{
572 MPQCData &data = *static_cast<MPQCData *>(_data);
573 // times obtain from key "mpqc" which should be the first
574 const int nregion = tim->nregion();
575 //std::cout << "There are " << nregion << " timed regions." << std::endl;
576 const char **region_names = new const char*[nregion];
577 tim->get_region_names(region_names);
578 // find "gather"
579 size_t gather_region = findTimerRegion(nregion, region_names, "gather");
580 size_t mpqc_region = findTimerRegion(nregion, region_names, "mpqc");
581 delete[] region_names;
582
583 // get timings
584 double *cpu_time = new double[nregion];
585 double *wall_time = new double[nregion];
586 double *flops = new double[nregion];
587 tim->get_cpu_times(cpu_time);
588 tim->get_wall_times(wall_time);
589 tim->get_flops(flops);
590 if (cpu_time != NULL) {
591 data.times.total_cputime = cpu_time[mpqc_region];
592 data.times.gather_cputime = cpu_time[gather_region];
593 }
594 if (wall_time != NULL) {
595 data.times.total_walltime = wall_time[mpqc_region];
596 data.times.gather_walltime = wall_time[gather_region];
597 }
598 if (flops != NULL) {
599 data.times.total_flops = flops[mpqc_region];
600 data.times.gather_flops = flops[gather_region];
601 }
602 delete[] cpu_time;
603 delete[] wall_time;
604 delete[] flops;
605}
606
Note: See TracBrowser for help on using the repository browser.