source: src/Jobs/InterfaceVMGJob.cpp@ 16b4fe

Action_Thermostats Add_AtomRandomPerturbation Add_FitFragmentPartialChargesAction Add_RotateAroundBondAction Add_SelectAtomByNameAction Added_ParseSaveFragmentResults AddingActions_SaveParseParticleParameters Adding_Graph_to_ChangeBondActions Adding_MD_integration_tests Adding_ParticleName_to_Atom Adding_StructOpt_integration_tests AtomFragments Automaking_mpqc_open AutomationFragmentation_failures Candidate_v1.5.4 Candidate_v1.6.0 Candidate_v1.6.1 Candidate_v1.7.0 ChangeBugEmailaddress ChangingTestPorts ChemicalSpaceEvaluator CombiningParticlePotentialParsing Combining_Subpackages Debian_Package_split Debian_package_split_molecuildergui_only Disabling_MemDebug Docu_Python_wait EmpiricalPotential_contain_HomologyGraph EmpiricalPotential_contain_HomologyGraph_documentation Enable_parallel_make_install Enhance_userguide Enhanced_StructuralOptimization Enhanced_StructuralOptimization_continued Example_ManyWaysToTranslateAtom Exclude_Hydrogens_annealWithBondGraph FitPartialCharges_GlobalError Fix_BoundInBox_CenterInBox_MoleculeActions Fix_ChargeSampling_PBC Fix_ChronosMutex Fix_FitPartialCharges Fix_FitPotential_needs_atomicnumbers Fix_ForceAnnealing Fix_IndependentFragmentGrids Fix_ParseParticles Fix_ParseParticles_split_forward_backward_Actions Fix_PopActions Fix_QtFragmentList_sorted_selection Fix_Restrictedkeyset_FragmentMolecule Fix_StatusMsg Fix_StepWorldTime_single_argument Fix_Verbose_Codepatterns Fix_fitting_potentials Fixes ForceAnnealing_goodresults ForceAnnealing_oldresults ForceAnnealing_tocheck ForceAnnealing_with_BondGraph ForceAnnealing_with_BondGraph_continued ForceAnnealing_with_BondGraph_continued_betteresults ForceAnnealing_with_BondGraph_contraction-expansion FragmentAction_writes_AtomFragments FragmentMolecule_checks_bonddegrees GeometryObjects Gui_Fixes Gui_displays_atomic_force_velocity ImplicitCharges IndependentFragmentGrids IndependentFragmentGrids_IndividualZeroInstances IndependentFragmentGrids_IntegrationTest IndependentFragmentGrids_Sole_NN_Calculation JobMarket_RobustOnKillsSegFaults JobMarket_StableWorkerPool JobMarket_unresolvable_hostname_fix MoreRobust_FragmentAutomation ODR_violation_mpqc_open PartialCharges_OrthogonalSummation PdbParser_setsAtomName PythonUI_with_named_parameters QtGui_reactivate_TimeChanged_changes Recreated_GuiChecks Rewrite_FitPartialCharges RotateToPrincipalAxisSystem_UndoRedo SaturateAtoms_findBestMatching SaturateAtoms_singleDegree StoppableMakroAction Subpackage_CodePatterns Subpackage_JobMarket Subpackage_LinearAlgebra Subpackage_levmar Subpackage_mpqc_open Subpackage_vmg Switchable_LogView ThirdParty_MPQC_rebuilt_buildsystem TrajectoryDependenant_MaxOrder TremoloParser_IncreasedPrecision TremoloParser_MultipleTimesteps TremoloParser_setsAtomName Ubuntu_1604_changes stable
Last change on this file since 16b4fe was 16b4fe, checked in by Frederik Heber <heber@…>, 13 years ago

DEBUG: Added writing grids as VTK files to InterfaceVMGJob.

  • writing both initial right-hand side and potential.
  • Property mode set to 100644
File size: 8.0 KB
Line 
1/*
2 * Project: MoleCuilder
3 * Description: creates and alters molecular systems
4 * Copyright (C) 2012 University of Bonn. All rights reserved.
5 *
6 *
7 * This file is part of MoleCuilder.
8 *
9 * MoleCuilder is free software: you can redistribute it and/or modify
10 * it under the terms of the GNU General Public License as published by
11 * the Free Software Foundation, either version 2 of the License, or
12 * (at your option) any later version.
13 *
14 * MoleCuilder is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 * GNU General Public License for more details.
18 *
19 * You should have received a copy of the GNU General Public License
20 * along with MoleCuilder. If not, see <http://www.gnu.org/licenses/>.
21 */
22
23/*
24 * InterfaceVMGJob.cpp
25 *
26 * Created on: 10.06.2012
27 * Author: Frederik Heber
28 */
29
30#ifdef HAVE_CONFIG_H
31#include <config.h>
32#endif
33
34#ifdef HAVE_MPI
35#include "mpi.h"
36#endif
37
38#include "CodePatterns/MemDebug.hpp"
39
40#include <cmath>
41#include <iostream>
42
43#include "CodePatterns/Log.hpp"
44
45#include "base/vector.hpp"
46#include "base/math.hpp"
47#ifdef HAVE_MPI
48#include "comm/comm_mpi.hpp"
49#else
50#include "comm/comm_serial.hpp"
51#endif
52#include "grid/grid.hpp"
53#include "grid/multigrid.hpp"
54
55#include "InterfaceVMGJob.hpp"
56
57using namespace VMG;
58using VMGInterfaces::InterfaceVMGJob;
59
60InterfaceVMGJob::InterfaceVMGJob(const std::vector< double > &_sampled_input,
61 VMGData &_returndata,
62 const std::vector< std::vector<double> > &_particle_positions,
63 const std::vector< double > &_particle_charges,
64 VMG::Boundary boundary,
65 int levelMin,
66 int levelMax,
67 const VMG::Vector &box_begin,
68 vmg_float box_end,
69 const int& near_field_cells,
70 int coarseningSteps,
71 double alpha) :
72 VMG::Interface(boundary, levelMin, levelMax,
73 box_begin, box_end, coarseningSteps, alpha),
74 spl(near_field_cells, Extent(MaxLevel()).MeshWidth().Max()),
75 sampled_input(_sampled_input),
76 returndata(_returndata),
77 level(levelMax)
78{
79 std::vector< std::vector<double> >::const_iterator positer = _particle_positions.begin();
80 std::vector<double>::const_iterator chargeiter = _particle_charges.begin();
81 double pos[3];
82 for (; positer != _particle_positions.end(); ++positer, ++chargeiter) {
83 ASSERT( (*positer).size() == 3,
84 "InterfaceVMGJob::InterfaceVMGJob() - particle "
85 +toString(distance(_particle_positions.begin(), positer))+" has not exactly 3 coordinates.");
86 for (size_t i=0;i<3;++i)
87 pos[i] = (*positer)[i];
88 particles.push_back(Particle::Particle(pos, *chargeiter));
89 }
90}
91
92InterfaceVMGJob::~InterfaceVMGJob()
93{}
94
95void InterfaceVMGJob::ImportRightHandSide(Multigrid& multigrid)
96{
97 Index i;
98 Vector pos;
99 // VMG::TempGrid *temp_grid = new VMG::TempGrid(129, 0, 0., 1.);
100
101 Grid& grid = multigrid(multigrid.MaxLevel());
102 grid.ClearBoundary();
103
104 /// 1. assign nuclei as smeared-out charges to the grid
105
106 /*
107 * Charge assignment on the grid
108 */
109#ifdef HAVE_MPI
110 CommMPI& comm = *dynamic_cast<CommMPI*>(VMG::MG::GetComm());
111#else
112 CommSerial& comm = *dynamic_cast<CommSerial*>(VMG::MG::GetComm());
113#endif
114 Grid& particle_grid = comm.GetParticleGrid();
115
116 particle_grid.Clear();
117
118 assert(particle_grid.Global().LocalSize().IsComponentwiseGreater(
119 VMG::MG::GetFactory().GetObjectStorageVal<int>("PARTICLE_NEAR_FIELD_CELLS")));
120
121 LOG(1, "INFO: Creating particle grid for " << particles.size() << " particles.");
122 for (std::list<Particle::Particle>::iterator iter = particles.begin();
123 iter != particles.end(); ++iter) {
124 LOG(2, "DEBUG: Current particle is at " << (*iter).Pos()
125 << " with charge " << (*iter).Charge() << ".");
126 spl.SetSpline(particle_grid, *iter);
127 }
128
129 double particle_charges = 0.0;
130 for (std::list<Particle::Particle>::iterator iter=particles.begin(); iter!=particles.end(); ++iter)
131 particle_charges += iter->Charge();
132 particle_charges = MG::GetComm()->GlobalSumRoot(particle_charges);
133 LOG(1, "INFO: Total charge is " << particle_charges << ".");
134
135 // Communicate charges over halo
136 comm.CommFromGhosts(particle_grid);
137
138 // Assign charge values to the right hand side
139 double normalization = 1.;
140 if (!particles.empty()) {
141 normalization =
142 particle_charges/(
143 pow(VMG::MG::GetFactory().GetObjectStorageVal<int>("PARTICLE_NEAR_FIELD_CELLS"), 3) *
144 (double)particles.size()
145 );
146 }
147 LOG(1, "INFO: Normalization is " << normalization << ".");
148 for (int i=0; i<grid.Local().Size().X(); ++i)
149 for (int j=0; j<grid.Local().Size().Y(); ++j)
150 for (int k=0; k<grid.Local().Size().Z(); ++k)
151 grid(grid.Local().Begin().X() + i,
152 grid.Local().Begin().Y() + j,
153 grid.Local().Begin().Z() + k) = normalization * 4.0 * VMG::Math::pi *
154 particle_grid.GetVal(particle_grid.Local().Begin().X() + i,
155 particle_grid.Local().Begin().Y() + j,
156 particle_grid.Local().Begin().Z() + k);
157
158 /// 2. add sampled electron density to the grid
159
160 const Index begin_local = grid.Global().LocalBegin() - grid.Local().HaloSize1();
161
162 LOG(1, "INFO: Mesh has extent " << grid.Extent().MeshWidth() << ".");
163 const int gridpoints = pow(2, level);
164 LOG(1, "INFO: gridpoints on finest level are " << gridpoints << ".");
165 assert( (grid.Extent().MeshWidth().X() * gridpoints) == 1 );
166 assert( (grid.Extent().MeshWidth().Y() * gridpoints) == 1 );
167 assert( (grid.Extent().MeshWidth().Z() * gridpoints) == 1 );
168 LOG(1, "INFO: "
169 << "X in [" << grid.Local().Begin().X() << "," << grid.Local().End().X() << "],"
170 << "Y in [" << grid.Local().Begin().Y() << "," << grid.Local().End().Y() << "],"
171 << "Z in [" << grid.Local().Begin().Z() << "," << grid.Local().End().Z() << "].");
172
173 const double element_volume =
174 grid.Extent().MeshWidth().X() * grid.Extent().MeshWidth().Y() * grid.Extent().MeshWidth().Z();
175 std::vector<double>::const_iterator sample_iter = sampled_input.begin();
176 for (i.X()=grid.Local().Begin().X(); i.X()<grid.Local().End().X(); ++i.X())
177 for (i.Y()=grid.Local().Begin().Y(); i.Y()<grid.Local().End().Y(); ++i.Y())
178 for (i.Z()=grid.Local().Begin().Z(); i.Z()<grid.Local().End().Z(); ++i.Z())
179 grid(i) -= (*sample_iter++) * element_volume; //temp_grid(i);
180 assert( sample_iter == sampled_input.end() );
181
182 Grid::iterator grid_iter;
183 double charge_sum = 0.0;
184 for (grid_iter=grid.Iterators().Local().Begin(); grid_iter!=grid.Iterators().Local().End(); ++grid_iter)
185 charge_sum += grid.GetVal(*grid_iter);
186// charge_sum = MG::GetComm()->GlobalSum(charge_sum);
187// comm.PrintStringOnce("Grid charge sum: %e", charge_sum);
188 std::cout << "Grid charge sum: " << charge_sum << std::endl;
189
190 // print input grid to vtk
191 VMG::MG::GetComm()->PrintGrid(grid, "grid after ImportRightHandSide");
192
193// delete temp_grid;
194}
195
196void InterfaceVMGJob::ExportSolution(Grid& grid)
197{
198 // grid now contains the sough-for potential
199 VMG::Comm* comm = VMG::MG::GetComm();
200
201 const Index begin_local = grid.Global().LocalBegin() - grid.Local().HaloSize1();
202 Index i;
203
204 // print output grid to vtk
205 comm->PrintGrid(grid, "grid after ExportSolution");
206
207 returndata.sampled_potential.sampled_grid.clear();
208 for (i.X()=grid.Local().Begin().X(); i.X()<grid.Local().End().X(); ++i.X())
209 for (i.Y()=grid.Local().Begin().Y(); i.Y()<grid.Local().End().Y(); ++i.Y())
210 for (i.Z()=grid.Local().Begin().Z(); i.Z()<grid.Local().End().Z(); ++i.Z()) {
211 returndata.sampled_potential.sampled_grid.push_back(grid(i));
212 }
213
214 Grid::iterator grid_iter;
215 double potential_sum = 0.0;
216 for (grid_iter=grid.Iterators().Local().Begin(); grid_iter!=grid.Iterators().Local().End(); ++grid_iter)
217 potential_sum += grid.GetVal(*grid_iter);
218// charge_sum = MG::GetComm()->GlobalSum(potential_sum);
219// comm.PrintStringOnce("Grid potential sum: %e", potential_sum);
220 std::cout << "Grid potential sum: " << potential_sum << std::endl;
221
222 //Particle::CommMPI& comm = *dynamic_cast<Particle::CommMPI*>(MG::GetComm());
223 returndata.e_long = VMG::MG::GetComm()->GlobalSumRoot(returndata.e_long);
224
225}
Note: See TracBrowser for help on using the repository browser.