source: src/Graph/BondGraph.cpp@ 440ac3

Action_Thermostats Add_AtomRandomPerturbation Add_RotateAroundBondAction Add_SelectAtomByNameAction Adding_Graph_to_ChangeBondActions Adding_MD_integration_tests Adding_StructOpt_integration_tests AutomationFragmentation_failures Candidate_v1.6.0 Candidate_v1.6.1 ChangeBugEmailaddress ChangingTestPorts ChemicalSpaceEvaluator Combining_Subpackages Debian_Package_split Debian_package_split_molecuildergui_only Disabling_MemDebug Docu_Python_wait EmpiricalPotential_contain_HomologyGraph_documentation Enable_parallel_make_install Enhance_userguide Enhanced_StructuralOptimization Enhanced_StructuralOptimization_continued Example_ManyWaysToTranslateAtom Exclude_Hydrogens_annealWithBondGraph FitPartialCharges_GlobalError Fix_ChronosMutex Fix_StatusMsg Fix_StepWorldTime_single_argument Fix_Verbose_Codepatterns ForceAnnealing_goodresults ForceAnnealing_oldresults ForceAnnealing_tocheck ForceAnnealing_with_BondGraph ForceAnnealing_with_BondGraph_continued ForceAnnealing_with_BondGraph_continued_betteresults ForceAnnealing_with_BondGraph_contraction-expansion GeometryObjects Gui_displays_atomic_force_velocity IndependentFragmentGrids_IntegrationTest JobMarket_RobustOnKillsSegFaults JobMarket_StableWorkerPool JobMarket_unresolvable_hostname_fix PartialCharges_OrthogonalSummation PythonUI_with_named_parameters QtGui_reactivate_TimeChanged_changes Recreated_GuiChecks RotateToPrincipalAxisSystem_UndoRedo StoppableMakroAction TremoloParser_IncreasedPrecision TremoloParser_MultipleTimesteps Ubuntu_1604_changes stable
Last change on this file since 440ac3 was 0763ce, checked in by Frederik Heber <heber@…>, 10 years ago

Added BondGraph::checkBondDegree, FragmentationAction only resets degrees when incorrect.

  • this fixes the bug where the molecular dynamics actions would flip the double bonds in an aromatic ring during the simulation steps because the bond degrees are reset even though the bond graph is present and should be re-used.
  • Property mode set to 100644
File size: 18.1 KB
Line 
1/*
2 * Project: MoleCuilder
3 * Description: creates and alters molecular systems
4 * Copyright (C) 2010-2012 University of Bonn. All rights reserved.
5 * Copyright (C) 2013 Frederik Heber. All rights reserved.
6 *
7 *
8 * This file is part of MoleCuilder.
9 *
10 * MoleCuilder is free software: you can redistribute it and/or modify
11 * it under the terms of the GNU General Public License as published by
12 * the Free Software Foundation, either version 2 of the License, or
13 * (at your option) any later version.
14 *
15 * MoleCuilder is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU General Public License for more details.
19 *
20 * You should have received a copy of the GNU General Public License
21 * along with MoleCuilder. If not, see <http://www.gnu.org/licenses/>.
22 */
23
24/*
25 * bondgraph.cpp
26 *
27 * Created on: Oct 29, 2009
28 * Author: heber
29 */
30
31// include config.h
32#ifdef HAVE_CONFIG_H
33#include <config.h>
34#endif
35
36// boost::graph uses placement new
37#include <boost/graph/adjacency_list.hpp>
38
39#include "CodePatterns/MemDebug.hpp"
40
41#include <algorithm>
42#include <boost/bimap.hpp>
43#include <boost/bind.hpp>
44#include <boost/foreach.hpp>
45#include <boost/function.hpp>
46#include <boost/graph/max_cardinality_matching.hpp>
47#include <deque>
48#include <iostream>
49
50#include "Atom/atom.hpp"
51#include "Bond/bond.hpp"
52#include "Graph/BondGraph.hpp"
53#include "Box.hpp"
54#include "Element/element.hpp"
55#include "CodePatterns/Info.hpp"
56#include "CodePatterns/Log.hpp"
57#include "CodePatterns/Range.hpp"
58#include "CodePatterns/Verbose.hpp"
59#include "molecule.hpp"
60#include "Element/periodentafel.hpp"
61#include "Fragmentation/MatrixContainer.hpp"
62#include "Graph/DepthFirstSearchAnalysis.hpp"
63#include "LinearAlgebra/Vector.hpp"
64#include "World.hpp"
65#include "WorldTime.hpp"
66
67const double BondGraph::BondThreshold = 0.4; //!< CSD threshold in bond check which is the width of the interval whose center is the sum of the covalent radii
68
69BondGraph::BondGraph() :
70 BondLengthMatrix(NULL),
71 IsAngstroem(true)
72{}
73
74BondGraph::BondGraph(bool IsA) :
75 BondLengthMatrix(NULL),
76 IsAngstroem(IsA)
77{}
78
79BondGraph::~BondGraph()
80{
81 CleanupBondLengthTable();
82}
83
84void BondGraph::CleanupBondLengthTable()
85{
86 if (BondLengthMatrix != NULL) {
87 delete(BondLengthMatrix);
88 }
89}
90
91bool BondGraph::LoadBondLengthTable(
92 std::istream &input)
93{
94 Info FunctionInfo(__func__);
95 bool status = true;
96 MatrixContainer *TempContainer = NULL;
97
98 // allocate MatrixContainer
99 if (BondLengthMatrix != NULL) {
100 LOG(1, "MatrixContainer for Bond length already present, removing.");
101 delete(BondLengthMatrix);
102 BondLengthMatrix = NULL;
103 }
104 TempContainer = new MatrixContainer;
105
106 // parse in matrix
107 if ((input.good()) && (status = TempContainer->ParseMatrix(input, 0, 1, 0))) {
108 LOG(1, "Parsing bond length matrix successful.");
109 } else {
110 ELOG(1, "Parsing bond length matrix failed.");
111 status = false;
112 }
113
114 if (status) // set to not NULL only if matrix was parsed
115 BondLengthMatrix = TempContainer;
116 else {
117 BondLengthMatrix = NULL;
118 delete(TempContainer);
119 }
120 return status;
121}
122
123double BondGraph::GetBondLength(
124 int firstZ,
125 int secondZ) const
126{
127 double return_length;
128 if ((firstZ < 0) || (firstZ >= (int)BondLengthMatrix->Matrix[0].size()))
129 return -1.;
130 if ((secondZ < 0) || (secondZ >= (int)BondLengthMatrix->Matrix[0][firstZ].size()))
131 return -1.;
132 if (BondLengthMatrix == NULL) {
133 return_length = -1.;
134 } else {
135 return_length = BondLengthMatrix->Matrix[0][firstZ][secondZ];
136 }
137 LOG(4, "INFO: Request for length between " << firstZ << " and "
138 << secondZ << ": " << return_length << ".");
139 return return_length;
140}
141
142range<double> BondGraph::CovalentMinMaxDistance(
143 const element * const Walker,
144 const element * const OtherWalker) const
145{
146 range<double> MinMaxDistance(0.,0.);
147 MinMaxDistance.first = OtherWalker->getCovalentRadius() + Walker->getCovalentRadius();
148 MinMaxDistance.first *= (IsAngstroem) ? 1. : 1. / AtomicLengthToAngstroem;
149 MinMaxDistance.last = MinMaxDistance.first + BondThreshold;
150 MinMaxDistance.first -= BondThreshold;
151 return MinMaxDistance;
152}
153
154range<double> BondGraph::BondLengthMatrixMinMaxDistance(
155 const element * const Walker,
156 const element * const OtherWalker) const
157{
158 ASSERT(BondLengthMatrix, "BondGraph::BondLengthMatrixMinMaxDistance() called without NULL BondLengthMatrix.");
159 ASSERT(Walker, "BondGraph::BondLengthMatrixMinMaxDistance() - illegal element given.");
160 ASSERT(OtherWalker, "BondGraph::BondLengthMatrixMinMaxDistance() - illegal other element given.");
161 range<double> MinMaxDistance(0.,0.);
162 MinMaxDistance.first = GetBondLength(Walker->getAtomicNumber()-1, OtherWalker->getAtomicNumber()-1);
163 MinMaxDistance.first *= (IsAngstroem) ? 1. : 1. / AtomicLengthToAngstroem;
164 MinMaxDistance.last = MinMaxDistance.first + BondThreshold;
165 MinMaxDistance.first -= BondThreshold;
166 return MinMaxDistance;
167}
168
169range<double> BondGraph::getMinMaxDistance(
170 const element * const Walker,
171 const element * const OtherWalker) const
172{
173 range<double> MinMaxDistance(0.,0.);
174 if (BondLengthMatrix == NULL) {// safety measure if no matrix has been parsed yet
175 LOG(5, "DEBUG: Using Covalent radii criterion for [min,max) distances.");
176 MinMaxDistance = CovalentMinMaxDistance(Walker, OtherWalker);
177 } else {
178 LOG(5, "DEBUG: Using tabulated bond table criterion for [min,max) distances.");
179 MinMaxDistance = BondLengthMatrixMinMaxDistance(Walker, OtherWalker);
180 }
181 return MinMaxDistance;
182}
183
184range<double> BondGraph::getMinMaxDistance(
185 const BondedParticle * const Walker,
186 const BondedParticle * const OtherWalker) const
187{
188 return getMinMaxDistance(Walker->getType(), OtherWalker->getType());
189}
190
191range<double> BondGraph::getMinMaxDistanceSquared(
192 const BondedParticle * const Walker,
193 const BondedParticle * const OtherWalker) const
194{
195 // use non-squared version
196 range<double> MinMaxDistance(getMinMaxDistance(Walker, OtherWalker));
197 // and square
198 MinMaxDistance.first *= MinMaxDistance.first;
199 MinMaxDistance.last *= MinMaxDistance.last;
200 return MinMaxDistance;
201}
202
203LinkedCell::LinkedCell_View BondGraph::getLinkedCell(const double max_distance) const
204{
205 return World::getInstance().getLinkedCell(max_distance);
206}
207
208std::set< const element *> BondGraph::getElementSetFromNumbers(const std::set<atomicNumber_t> &Set) const
209{
210 std::set< const element *> PresentElements;
211 for(std::set<atomicNumber_t>::const_iterator iter = Set.begin(); iter != Set.end(); ++iter)
212 PresentElements.insert( World::getInstance().getPeriode()->FindElement(*iter) );
213 return PresentElements;
214}
215
216Box &BondGraph::getDomain() const
217{
218 return World::getInstance().getDomain();
219}
220
221unsigned int BondGraph::getTime() const
222{
223 return WorldTime::getTime();
224}
225
226bool BondGraph::operator==(const BondGraph &other) const
227{
228 if (IsAngstroem != other.IsAngstroem)
229 return false;
230 if (((BondLengthMatrix != NULL) && (other.BondLengthMatrix == NULL))
231 || ((BondLengthMatrix == NULL) && (other.BondLengthMatrix != NULL)))
232 return false;
233 if ((BondLengthMatrix != NULL) && (other.BondLengthMatrix != NULL))
234 if (*BondLengthMatrix != *other.BondLengthMatrix)
235 return false;
236 return true;
237}
238
239/** Corrects the bond degree by one at most if necessary.
240 *
241 * We are constraint to bonds in \a egdes, if the candidate bond is in edges, we
242 * just don't increment its bond degree. We do not choose an alternative for this
243 * atom.
244 *
245 * \param edges only these edges can be updated
246 * \return number of corrections done
247 */
248int BondGraph::CorrectBondDegree(atom *_atom, const std::set<bond::ptr>& edges) const
249{
250 int NoBonds = 0;
251 int OtherNoBonds = 0;
252 int FalseBondDegree = 0;
253 int CandidateDeltaNoBonds = -1;
254 atom *OtherWalker = NULL;
255 bond::ptr CandidateBond;
256
257 NoBonds = _atom->CountBonds();
258 LOG(3, "Walker " << *_atom << ": " << (int)_atom->getType()->getNoValenceOrbitals() << " > " << NoBonds << "?");
259 const int DeltaNoBonds = (int)(_atom->getType()->getNoValenceOrbitals()) - NoBonds;
260 if (DeltaNoBonds > 0) { // we have a mismatch, check all bonding partners for mismatch
261 const BondList& ListOfBonds = _atom->getListOfBonds();
262 for (BondList::const_iterator Runner = ListOfBonds.begin(); Runner != ListOfBonds.end(); (++Runner)) {
263 OtherWalker = (*Runner)->GetOtherAtom(_atom);
264 OtherNoBonds = OtherWalker->CountBonds();
265 const int OtherDeltaNoBonds = (int)(OtherWalker->getType()->getNoValenceOrbitals()) - OtherNoBonds;
266 LOG(3, "OtherWalker " << *OtherWalker << ": " << (int)OtherWalker->getType()->getNoValenceOrbitals() << " > " << OtherNoBonds << "?");
267 if (OtherDeltaNoBonds > 0) { // check if possible candidate
268 const BondList& OtherListOfBonds = OtherWalker->getListOfBonds();
269 if ((CandidateBond == NULL)
270 || (ListOfBonds.size() > OtherListOfBonds.size()) // pick the one with fewer number of bonds first
271 || ((CandidateDeltaNoBonds < 0) || (CandidateDeltaNoBonds > OtherDeltaNoBonds)) ) // pick the one closer to fulfilling its valence first
272 {
273 CandidateDeltaNoBonds = OtherDeltaNoBonds;
274 CandidateBond = (*Runner);
275 LOG(3, "New candidate is " << *CandidateBond << ".");
276 }
277 }
278 }
279 if ((CandidateBond != NULL)) {
280 if (edges.find(CandidateBond) != edges.end()) {
281 CandidateBond->setDegree(CandidateBond->getDegree()+1);
282 LOG(2, "Increased bond degree for bond " << *CandidateBond << ".");
283 } else
284 LOG(2, "Bond " << *CandidateBond << " is not in edges.");
285 } else {
286 ELOG(2, "Could not find correct degree for atom " << *_atom << ".");
287 FalseBondDegree++;
288 }
289 }
290 return FalseBondDegree;
291}
292
293/** Sets the weight of all connected bonds to one.
294 */
295void BondGraph::resetBondDegree(atom *_atom) const
296{
297 const BondList &ListOfBonds = _atom->getListOfBonds();
298 for (BondList::const_iterator BondRunner = ListOfBonds.begin();
299 BondRunner != ListOfBonds.end();
300 ++BondRunner)
301 (*BondRunner)->setDegree(1);
302}
303
304std::set<bond::ptr> BondGraph::getBackEdges() const
305{
306 DepthFirstSearchAnalysis DFS;
307 DFS();
308
309 const std::deque<bond::ptr>& backedgestack = DFS.getBackEdgeStack();
310 std::set<bond::ptr> backedges(backedgestack.begin(), backedgestack.end());
311
312 return backedges;
313}
314
315std::set<bond::ptr> BondGraph::getTreeEdges() const
316{
317 const std::set<bond::ptr> backedges = getBackEdges();
318 std::set<bond::ptr> alledges;
319 const World& world = World::getInstance();
320 for(World::AtomConstIterator iter = world.getAtomIter();
321 iter != world.atomEnd(); ++iter) {
322 const BondList &ListOfBonds = (*iter)->getListOfBonds();
323 alledges.insert(ListOfBonds.begin(), ListOfBonds.end());
324 }
325 std::set<bond::ptr> treeedges;
326 std::set_difference(
327 alledges.begin(), alledges.end(),
328 backedges.begin(), backedges.end(),
329 std::inserter(treeedges, treeedges.begin()));
330 return treeedges;
331}
332
333struct hasDelta {
334 bool operator()(atom *_atom) {
335 const double delta =
336 _atom->getType()->getNoValenceOrbitals() - _atom->CountBonds();
337 return (fabs(delta) > 0);
338 }
339};
340
341typedef std::set<atom *> deltaatoms_t;
342typedef std::set<bond::ptr> deltaedges_t;
343
344static int gatherAllDeltaAtoms(
345 const deltaatoms_t &allatoms,
346 deltaatoms_t &deltaatoms)
347{
348 static hasDelta delta;
349 deltaatoms.clear();
350 for (deltaatoms_t::const_iterator iter = allatoms.begin();
351 iter != allatoms.end(); ++iter) {
352 if (delta(*iter))
353 deltaatoms.insert(*iter);
354 }
355 return deltaatoms.size();
356}
357
358typedef boost::bimap<int,atom*> AtomIndexLookup_t;
359
360static AtomIndexLookup_t createAtomIndexLookup(
361 const deltaedges_t &edges)
362{
363 AtomIndexLookup_t AtomIndexLookup;
364 size_t index = 0;
365 for (deltaedges_t::const_iterator edgeiter = edges.begin();
366 edgeiter != edges.end(); ++edgeiter) {
367 const bond::ptr & _bond = *edgeiter;
368
369 // insert left into lookup
370 std::pair< AtomIndexLookup_t::iterator, bool > inserter =
371 AtomIndexLookup.insert( AtomIndexLookup_t::value_type(index, _bond->leftatom) );
372 if (inserter.second)
373 ++index;
374
375 // insert right into lookup
376 inserter = AtomIndexLookup.insert( AtomIndexLookup_t::value_type(index, _bond->rightatom) );
377 if (inserter.second)
378 ++index;
379 }
380 return AtomIndexLookup;
381}
382
383typedef std::map< std::pair<atom *, atom*>, bond::ptr> BondLookup_t;
384static BondLookup_t createBondLookup(
385 const deltaedges_t &edges)
386{
387 BondLookup_t BondLookup;
388 for (deltaedges_t::const_iterator edgeiter = edges.begin();
389 edgeiter != edges.end(); ++edgeiter) {
390 const bond::ptr & _bond = *edgeiter;
391
392 // insert bond into pair lookup
393 BondLookup.insert(
394 std::make_pair(
395 std::make_pair(_bond->leftatom, _bond->rightatom),
396 _bond)
397 );
398 }
399 return BondLookup;
400}
401
402typedef boost::adjacency_list<boost::vecS, boost::vecS, boost::undirectedS> graph_t;
403typedef std::vector<boost::graph_traits<graph_t>::vertex_descriptor> mate_t;
404
405static void fillEdgesIntoMatching(
406 graph_t &g,
407 mate_t &mate,
408 const AtomIndexLookup_t &AtomIndexLookup,
409 const BondLookup_t &BondLookup,
410 deltaedges_t &matching
411 )
412{
413 matching.clear();
414 boost::graph_traits<graph_t>::vertex_iterator vi, vi_end;
415 for(tie(vi,vi_end) = boost::vertices(g); vi != vi_end; ++vi)
416 if (mate[*vi] != boost::graph_traits<graph_t>::null_vertex() && *vi < mate[*vi]) {
417 atom * leftatom = AtomIndexLookup.left.at(*vi);
418 atom * rightatom = AtomIndexLookup.left.at(mate[*vi]);
419 std::pair<atom *,atom *> AtomPair(leftatom, rightatom);
420 std::pair<atom *,atom *> reverseAtomPair(rightatom, leftatom);
421 BondLookup_t::const_iterator iter = BondLookup.find(AtomPair);
422 if (iter != BondLookup.end()) {
423 matching.insert(iter->second);
424 } else {
425 iter = BondLookup.find(reverseAtomPair);
426 if (iter != BondLookup.end()) {
427 matching.insert(iter->second);
428 } else
429 ASSERT(0, "fillEdgesIntoMatching() - cannot find ("+toString(*vi)+","
430 +toString(mate[*vi])+") in BondLookup.");
431 }
432 }
433}
434
435static bool createMatching(
436 deltaedges_t &deltaedges,
437 deltaedges_t &matching
438 )
439{
440 // gather all vertices for graph in a lookup structure
441 // to translate the graphs indices to atoms and also to associated bonds
442 AtomIndexLookup_t AtomIndexLookup = createAtomIndexLookup(deltaedges);
443 BondLookup_t BondLookup = createBondLookup(deltaedges);
444 const int n_vertices = AtomIndexLookup.size();
445
446 // construct graph
447 graph_t g(n_vertices);
448
449 // add edges to graph
450 for (deltaedges_t::const_iterator edgeiter = deltaedges.begin();
451 edgeiter != deltaedges.end(); ++edgeiter) {
452 const bond::ptr & _bond = *edgeiter;
453 boost::add_edge(
454 AtomIndexLookup.right.at(_bond->leftatom),
455 AtomIndexLookup.right.at(_bond->rightatom),
456 g);
457 }
458
459 // mate structure contains unique edge partner to every vertex in matching
460 mate_t mate(n_vertices);
461
462 // get maximum matching over given edges
463 bool success = boost::checked_edmonds_maximum_cardinality_matching(g, &mate[0]);
464
465 if (success) {
466 LOG(1, "STATUS: Found a matching of size " << matching_size(g, &mate[0]) << ".");
467 fillEdgesIntoMatching(g, mate, AtomIndexLookup, BondLookup, matching);
468 } else {
469 ELOG(2, "Failed to find maximum matching.");
470 }
471
472 return success;
473}
474
475bool BondGraph::checkBondDegree(const deltaatoms_t &allatoms) const
476{
477 deltaatoms_t deltaatoms;
478 return (gatherAllDeltaAtoms(allatoms, deltaatoms) == 0);
479}
480
481
482int BondGraph::calculateBondDegree(const deltaatoms_t &allatoms) const
483{
484 deltaatoms_t deltaatoms;
485 deltaedges_t deltaedges;
486 deltaedges_t matching;
487
488 LOG(1, "STATUS: Calculating bond degrees.");
489 if (DoLog(2)) {
490 std::stringstream output;
491 output << "INFO: All atoms are: ";
492 BOOST_FOREACH( atom *_atom, allatoms) {
493 output << *_atom << " ";
494 }
495 LOG(2, output.str());
496 }
497
498 size_t IterCount = 0;
499 // 1. gather all atoms with delta > 0
500 while ((gatherAllDeltaAtoms(allatoms, deltaatoms) != 0) && (IterCount < 100)) {
501 if (DoLog(3)) {
502 std::stringstream output;
503 output << "DEBUG: Current delta atoms are: ";
504 BOOST_FOREACH( atom *_atom, deltaatoms) {
505 output << *_atom << " ";
506 }
507 LOG(3, output.str());
508 }
509
510 // 2. gather all edges that connect atoms with both(!) having delta > 0
511 deltaedges.clear();
512 for (deltaatoms_t::const_iterator atomiter = deltaatoms.begin();
513 atomiter != deltaatoms.end(); ++atomiter) {
514 const atom * const _atom = *atomiter;
515 const BondList& ListOfBonds = (*atomiter)->getListOfBonds();
516 for (BondList::const_iterator bonditer = ListOfBonds.begin();
517 bonditer != ListOfBonds.end();++bonditer) {
518 const bond::ptr _bond = *bonditer;
519 if ((_bond->leftatom == _atom) && (deltaatoms.find(_bond->rightatom) != deltaatoms.end()))
520 deltaedges.insert(_bond);
521 else if ((_bond->rightatom == _atom) && (deltaatoms.find(_bond->leftatom) != deltaatoms.end()))
522 deltaedges.insert(_bond);
523 }
524 }
525 if (DoLog(3)) {
526 std::stringstream output;
527 output << "DEBUG: Current delta bonds are: ";
528 BOOST_FOREACH( bond::ptr _bond, deltaedges) {
529 output << *_bond << " ";
530 }
531 LOG(3, output.str());
532 }
533 if (deltaedges.empty())
534 break;
535
536 // 3. build matching over these edges
537 if (!createMatching(deltaedges, matching))
538 break;
539 if (DoLog(2)) {
540 std::stringstream output;
541 output << "DEBUG: Resulting matching is: ";
542 BOOST_FOREACH( bond::ptr _bond, matching) {
543 output << *_bond << " ";
544 }
545 LOG(2, output.str());
546 }
547
548 // 4. increase degree for each edge of the matching
549 LOG(2, "DEBUG: Increasing bond degrees by one.");
550 for (deltaedges_t::iterator edgeiter = matching.begin();
551 edgeiter != matching.end(); ++edgeiter) {
552 (*edgeiter)->setDegree( (*edgeiter)->getDegree()+1 );
553 }
554
555 // 6. stop after a given number of iterations or when all atoms are happy.
556 ++IterCount;
557 };
558
559 // check whether we still have deltaatoms
560 {
561 hasDelta delta;
562 for (deltaatoms_t::const_iterator iter = allatoms.begin();
563 iter != allatoms.end(); ++iter)
564 if (delta(*iter))
565 ELOG(2, "Could not find satisfy charge neutrality for atom " << **iter << ".");
566 }
567
568 return deltaedges.size();
569}
Note: See TracBrowser for help on using the repository browser.