ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-4/src/brains/SimCreator.cpp
Revision: 1807
Committed: Tue Nov 30 22:43:51 2004 UTC (19 years, 9 months ago) by tim
File size: 18977 byte(s)
Log Message:
brains get built

File Contents

# User Rev Content
1 tim 1703 /*
2     * Copyright (C) 2000-2004 Object Oriented Parallel Simulation Engine (OOPSE) project
3     *
4     * Contact: oopse@oopse.org
5     *
6     * This program is free software; you can redistribute it and/or
7     * modify it under the terms of the GNU Lesser General Public License
8     * as published by the Free Software Foundation; either version 2.1
9     * of the License, or (at your option) any later version.
10     * All we ask is that proper credit is given for our work, which includes
11     * - but is not limited to - adding the above copyright notice to the beginning
12     * of your source code files, and to any copyright notice that you may distribute
13     * with programs based on this work.
14     *
15     * This program 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 Lesser General Public License for more details.
19     *
20     * You should have received a copy of the GNU Lesser General Public License
21     * along with this program; if not, write to the Free Software
22     * Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
23     *
24     */
25    
26 tim 1725 /**
27     * @file SimCreator.cpp
28     * @author tlin
29     * @date 11/03/2004
30     * @time 13:51am
31     * @version 1.0
32     */
33    
34 tim 1807 #include <sprng.h>
35    
36     #include "brains/MoleculeCreator.hpp"
37 tim 1725 #include "brains/SimCreator.hpp"
38 tim 1807 #include "brains/SimSnapshotManager.hpp"
39     #include "io/DumpReader.hpp"
40     #include "io/parse_me.h"
41     #include "UseTheForce/ForceFieldFactory.hpp"
42     #include "utils/simError.h"
43     #include "utils/StringUtils.hpp"
44     #ifdef IS_MPI
45     #include "io/mpiBASS.h"
46     #endif
47 tim 1733
48 tim 1703 namespace oopse {
49    
50 tim 1807 void SimCreator::parseFile(const std::string mdFileName, MakeStamps* stamps, Globals* globals){
51 tim 1733
52 tim 1712 #ifdef IS_MPI
53 tim 1733
54     if (worldRank == 0) {
55 tim 1712 #endif // is_mpi
56    
57 tim 1733 globals->initalize();
58     set_interface_stamps(stamps, globals);
59 tim 1712
60 tim 1733 #ifdef IS_MPI
61 tim 1712
62 tim 1733 mpiEventInit();
63    
64 tim 1712 #endif
65    
66 tim 1733 yacc_BASS(mdFileName.c_str());
67 tim 1712
68     #ifdef IS_MPI
69 tim 1733
70     throwMPIEvent(NULL);
71     } else {
72     set_interface_stamps(stamps, globals);
73     mpiEventInit();
74     MPIcheckPoint();
75     mpiEventLoop();
76     }
77    
78 tim 1712 #endif
79    
80 tim 1703 }
81    
82 tim 1733 SimInfo* SimCreator::createSim(const std::string & mdFileName) {
83 tim 1712
84 tim 1733 MakeStamps * stamps = new MakeStamps();
85 tim 1712
86 tim 1733 Globals * globals = new Globals();
87    
88 tim 1712 //parse meta-data file
89 tim 1807 parseFile(mdFileName, stamps, globals);
90 tim 1712
91     //create the force field
92 tim 1807 ForceField * ff = ForceFieldFactory::getInstance()->createForceField(
93 tim 1733 globals->getForceField());
94 tim 1740
95 tim 1733 if (ff == NULL) {
96 tim 1807 sprintf(painCave.errMsg, "ForceField Factory can not create %s force field\n",
97 tim 1733 globals->getForceField());
98     painCave.isFatal = 1;
99     simError();
100     }
101    
102 tim 1807 std::string forcefieldFileName;
103     forcefieldFileName = ff->getForceFieldFileName();
104    
105     if (globals->haveForceFieldVariant()) {
106     //If the force field has variant, the variant force field name will be
107     //Base.variant.frc. For exampel EAM.u6.frc
108    
109     std::string variant = globals->getForceFieldVariant();
110    
111     std::string::size_type pos = forcefieldFileName.rfind(".frc");
112     variant = "." + variant;
113     if (pos != std::string::npos) {
114     forcefieldFileName.insert(pos, variant);
115     } else {
116     //If the default force field file name does not containt .frc suffix, just append the .variant
117     forcefieldFileName.append(variant);
118     }
119     }
120    
121     ff->parse(forcefieldFileName);
122    
123 tim 1733 //extract the molecule stamps
124     std::vector < std::pair<MoleculeStamp *, int> > moleculeStampPairs;
125     compList(stamps, globals, moleculeStampPairs);
126    
127 tim 1719 //create SimInfo
128 tim 1733 SimInfo * info = new SimInfo(moleculeStampPairs, ff, globals);
129 tim 1719
130 tim 1733 //gather parameters (SimCreator only retrieves part of the parameters)
131 tim 1807 gatherParameters(info, mdFileName);
132 tim 1733
133 tim 1712 //divide the molecules and determine the global index of molecules
134 tim 1807 #ifdef IS_MPI
135 tim 1733 divideMolecules(info);
136 tim 1807 #endif
137 tim 1733
138 tim 1712 //create the molecules
139 tim 1719 createMolecules(info);
140 tim 1712
141 tim 1807
142     //allocate memory for DataStorage(circular reference, need to break it)
143     info->setSnapshotManager(new SimSnapshotManager(info));
144    
145 tim 1725 //set the global index of atoms, rigidbodies and cutoffgroups (only need to be set once, the
146     //global index will never change again). Local indices of atoms and rigidbodies are already set by
147 tim 1733 //MoleculeCreator class which actually delegates the responsibility to LocalIndexManager.
148 tim 1807 setGlobalIndex(info);
149 tim 1733
150 tim 1807
151 tim 1733
152 tim 1725 //load initial coordinates, some extra information are pushed into SimInfo's property map ( such as
153     //eta, chi for NPT integrator)
154 tim 1733 loadCoordinates(info);
155 tim 1719 return info;
156 tim 1733 }
157 tim 1712
158 tim 1807 void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) {
159 tim 1725
160     //setup seed for random number generator
161     int seedValue;
162 tim 1733 Globals * globals = info->getGlobals();
163    
164     if (globals->haveSeed()) {
165 tim 1725 seedValue = globals->getSeed();
166    
167 tim 1733 if (seedValue / 1000000000 == 0) {
168 tim 1725 sprintf(painCave.errMsg,
169 tim 1733 "Seed for sprng library should contain at least 9 digits\n"
170     "OOPSE will generate a seed for user\n");
171    
172 tim 1725 painCave.isFatal = 0;
173     simError();
174    
175 tim 1733 //using seed generated by system instead of invalid seed set by user
176    
177 tim 1725 #ifndef IS_MPI
178 tim 1733
179 tim 1725 seedValue = make_sprng_seed();
180 tim 1733
181     #else
182    
183     if (worldRank == 0) {
184 tim 1725 seedValue = make_sprng_seed();
185     }
186 tim 1733
187     MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
188    
189     #endif
190    
191 tim 1725 } //end if (seedValue /1000000000 == 0)
192 tim 1733 } else {
193    
194 tim 1725 #ifndef IS_MPI
195 tim 1733
196 tim 1725 seedValue = make_sprng_seed();
197 tim 1733
198     #else
199    
200     if (worldRank == 0) {
201 tim 1725 seedValue = make_sprng_seed();
202     }
203 tim 1733
204     MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
205    
206 tim 1725 #endif
207 tim 1733
208     } //end of globals->haveSeed()
209    
210 tim 1725 info->setSeed(seedValue);
211    
212 tim 1733
213     //figure out the ouput file names
214 tim 1725 std::string prefix;
215 tim 1733
216 tim 1725 #ifdef IS_MPI
217 tim 1733
218     if (worldRank == 0) {
219 tim 1725 #endif // is_mpi
220 tim 1733
221     if (globals->haveFinalConfig()) {
222 tim 1807 prefix = getPrefix(globals->getFinalConfig());
223 tim 1725 } else {
224 tim 1807 prefix = getPrefix(mdfile);
225 tim 1725 }
226 tim 1733
227     info->setFinalConfigFileName(prefix + ".eor");
228 tim 1725 info->setDumpFileName(prefix + ".dump");
229     info->setStatFileName(prefix + ".stat");
230    
231     #ifdef IS_MPI
232 tim 1733
233 tim 1725 }
234    
235 tim 1733 #endif
236 tim 1725
237 tim 1712 }
238    
239     #ifdef IS_MPI
240 tim 1733 void SimCreator::divideMolecules(SimInfo *info) {
241 tim 1727 double numerator;
242     double denominator;
243     double precast;
244     double x;
245     double y;
246     double a;
247     int old_atoms;
248     int add_atoms;
249     int new_atoms;
250     int nTarget;
251     int done;
252     int i;
253     int j;
254     int loops;
255     int which_proc;
256 tim 1733 int nProcessors;
257     int nGlobalMols;
258     std::vector < int > atomsPerProc;
259     randomSPRNG myRandom(info->getSeed());
260     int * molToProcMap = info->getMolToProcMapPointer();
261 tim 1727
262 tim 1733 MPI_Comm_size(MPI_COMM_WORLD, &nProcessors);
263     nGlobalMols = info->getNGlobalMolecules();
264 tim 1727
265 tim 1733 if (nProcessors > nGlobalMols) {
266 tim 1727 sprintf(painCave.errMsg,
267     "nProcessors (%d) > nMol (%d)\n"
268     "\tThe number of processors is larger than\n"
269     "\tthe number of molecules. This will not result in a \n"
270     "\tusable division of atoms for force decomposition.\n"
271     "\tEither try a smaller number of processors, or run the\n"
272 tim 1733 "\tsingle-processor version of OOPSE.\n", nProcessors, nGlobalMols);
273 tim 1727
274     painCave.isFatal = 1;
275     simError();
276     }
277    
278 tim 1733 a = 3.0 * nGlobalMols / info->getNGlobalAtoms();
279 tim 1727
280 tim 1733 //initialize atomsPerProc
281     atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0);
282 tim 1727
283 tim 1733 for(i = 0; i < nGlobalMols; i++) {
284 tim 1727 // default to an error condition:
285 tim 1733 molToProcMap[i] = -1;
286 tim 1727 }
287    
288 tim 1733 if (worldRank == 0) {
289     numerator = info->getNGlobalAtoms();
290     denominator = nProcessors;
291 tim 1727 precast = numerator / denominator;
292     nTarget = (int)(precast + 0.5);
293    
294 tim 1733 for(i = 0; i < nGlobalMols; i++) {
295 tim 1727 done = 0;
296     loops = 0;
297    
298     while (!done) {
299     loops++;
300    
301     // Pick a processor at random
302    
303 tim 1733 which_proc = (int) (myRandom.getRandom() * nProcessors);
304 tim 1727
305 tim 1733 //get the molecule stamp first
306     int stampId = info->getMoleculeStampId(i);
307     MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId);
308 tim 1727
309 tim 1733 // How many atoms does this processor have so far?
310     old_atoms = atomsPerProc[which_proc];
311     add_atoms = moleculeStamp->getNAtoms();
312 tim 1727 new_atoms = old_atoms + add_atoms;
313    
314     // If we've been through this loop too many times, we need
315     // to just give up and assign the molecule to this processor
316     // and be done with it.
317    
318     if (loops > 100) {
319 tim 1733 sprintf(painCave.errMsg,
320     "I've tried 100 times to assign molecule %d to a "
321     " processor, but can't find a good spot.\n"
322     "I'm assigning it at random to processor %d.\n",
323     i, which_proc);
324 tim 1727
325     painCave.isFatal = 0;
326     simError();
327    
328 tim 1733 molToProcMap[i] = which_proc;
329     atomsPerProc[which_proc] += add_atoms;
330 tim 1727
331     done = 1;
332     continue;
333     }
334    
335     // If we can add this molecule to this processor without sending
336     // it above nTarget, then go ahead and do it:
337    
338     if (new_atoms <= nTarget) {
339 tim 1733 molToProcMap[i] = which_proc;
340     atomsPerProc[which_proc] += add_atoms;
341 tim 1727
342     done = 1;
343     continue;
344     }
345    
346     // The only situation left is when new_atoms > nTarget. We
347     // want to accept this with some probability that dies off the
348     // farther we are from nTarget
349    
350     // roughly: x = new_atoms - nTarget
351     // Pacc(x) = exp(- a * x)
352     // where a = penalty / (average atoms per molecule)
353    
354     x = (double)(new_atoms - nTarget);
355 tim 1733 y = myRandom.getRandom();
356 tim 1727
357     if (y < exp(- a * x)) {
358 tim 1733 molToProcMap[i] = which_proc;
359     atomsPerProc[which_proc] += add_atoms;
360 tim 1727
361     done = 1;
362     continue;
363 tim 1733 } else {
364     continue;
365     }
366 tim 1727 }
367     }
368    
369     // Spray out this nonsense to all other processors:
370    
371 tim 1733 MPI_Bcast(molToProcMap, nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
372 tim 1727 } else {
373    
374 tim 1733 // Listen to your marching orders from processor 0:
375 tim 1727
376 tim 1733 MPI_Bcast(molToProcMap, nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
377 tim 1727 }
378    
379     sprintf(checkPointMsg,
380     "Successfully divided the molecules among the processors.\n");
381     MPIcheckPoint();
382 tim 1712 }
383 tim 1733
384 tim 1713 #endif
385 tim 1712
386 tim 1733 void SimCreator::createMolecules(SimInfo *info) {
387 tim 1725 MoleculeCreator molCreator;
388     int stampId;
389 tim 1713
390 tim 1733 for(int i = 0; i < info->getNGlobalMolecules(); i++) {
391 tim 1713
392 tim 1733 #ifdef IS_MPI
393    
394     if (info->getMolToProc(i) == worldRank) {
395     #endif
396    
397 tim 1725 stampId = info->getMoleculeStampId(i);
398 tim 1807 Molecule * mol = molCreator.createMolecule(info->getForceField(), info->getMoleculeStamp(stampId),
399     stampId, i, info->getLocalIndexManager());
400 tim 1733
401 tim 1719 info->addMolecule(mol);
402 tim 1733
403     #ifdef IS_MPI
404    
405 tim 1719 }
406 tim 1733
407     #endif
408    
409     } //end for(int i=0)
410 tim 1713 }
411    
412 tim 1775 void SimCreator::compList(MakeStamps *stamps, Globals* globals,
413 tim 1733 std::vector < std::pair<MoleculeStamp *, int> > &moleculeStampPairs) {
414 tim 1719 int i;
415 tim 1733 char * id;
416     MoleculeStamp * currentStamp;
417 tim 1807 Component** the_components = globals->getComponents();
418 tim 1733 int n_components = globals->getNComponents();
419 tim 1719
420 tim 1733 if (!globals->haveNMol()) {
421 tim 1719 // we don't have the total number of molecules, so we assume it is
422     // given in each component
423    
424 tim 1733 for(i = 0; i < n_components; i++) {
425     if (!the_components[i]->haveNMol()) {
426 tim 1719 // we have a problem
427     sprintf(painCave.errMsg,
428 tim 1807 "SimCreator Error. No global NMol or component NMol given.\n"
429 tim 1733 "\tCannot calculate the number of atoms.\n");
430    
431 tim 1719 painCave.isFatal = 1;
432     simError();
433     }
434    
435     id = the_components[i]->getType();
436 tim 1807 currentStamp = (stamps->extractMolStamp(id))->getStamp();
437 tim 1733
438     if (currentStamp == NULL) {
439 tim 1719 sprintf(painCave.errMsg,
440 tim 1807 "SimCreator error: Component \"%s\" was not found in the "
441 tim 1733 "list of declared molecules\n", id);
442    
443 tim 1719 painCave.isFatal = 1;
444     simError();
445 tim 1733 }
446 tim 1719
447 tim 1733 moleculeStampPairs.push_back(
448 tim 1807 make_pair(currentStamp, the_components[i]->getNMol()));
449 tim 1719 } //end for (i = 0; i < n_components; i++)
450 tim 1733 } else {
451     sprintf(painCave.errMsg, "SimSetup error.\n"
452     "\tSorry, the ability to specify total"
453     " nMols and then give molfractions in the components\n"
454     "\tis not currently supported."
455     " Please give nMol in the components.\n");
456    
457 tim 1719 painCave.isFatal = 1;
458     simError();
459     }
460 tim 1733
461 tim 1719 #ifdef IS_MPI
462 tim 1733
463     strcpy(checkPointMsg, "Component stamps successfully extracted\n");
464     MPIcheckPoint();
465    
466 tim 1719 #endif // is_mpi
467 tim 1733
468 tim 1719 }
469    
470 tim 1733 void SimCreator::setGlobalIndex(SimInfo *info) {
471 tim 1807 SimInfo::MoleculeIterator mi;
472     Molecule::AtomIterator ai;
473     Molecule::RigidBodyIterator ri;
474     Molecule::CutoffGroupIterator ci;
475 tim 1733 Molecule * mol;
476     Atom * atom;
477     RigidBody * rb;
478     CutoffGroup * cg;
479 tim 1725 int beginAtomIndex;
480     int beginRigidBodyIndex;
481     int beginCutoffGroupIndex;
482 tim 1719
483 tim 1733 #ifndef IS_MPI
484    
485 tim 1725 beginAtomIndex = 0;
486     beginRigidBodyIndex = 0;
487     beginCutoffGroupIndex = 0;
488 tim 1733
489 tim 1725 #else
490 tim 1733
491 tim 1725 int nproc;
492     int myNode;
493    
494     myNode = worldRank;
495     MPI_Comm_size(MPI_COMM_WORLD, &nproc);
496    
497 tim 1733 std::vector < int > tmpAtomsInProc(nproc, 0);
498     std::vector < int > tmpRigidBodiesInProc(nproc, 0);
499     std::vector < int > tmpCutoffGroupsInProc(nproc, 0);
500     std::vector < int > NumAtomsInProc(nproc, 0);
501     std::vector < int > NumRigidBodiesInProc(nproc, 0);
502     std::vector < int > NumCutoffGroupsInProc(nproc, 0);
503    
504     tmpAtomsInProc[myNode] = info->getNAtoms();
505     tmpRigidBodiesInProc[myNode] = info->getNRigidBodies();
506     tmpCutoffGroupsInProc[myNode] = info->getNCutoffGroups();
507    
508 tim 1725 //do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups
509 tim 1733 MPI_Allreduce(&tmpAtomsInProc[0], &NumAtomsInProc[0], nproc, MPI_INT,
510     MPI_SUM, MPI_COMM_WORLD);
511     MPI_Allreduce(&tmpRigidBodiesInProc[0], &NumRigidBodiesInProc[0], nproc,
512     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
513     MPI_Allreduce(&tmpCutoffGroupsInProc[0], &NumCutoffGroupsInProc[0], nproc,
514     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
515 tim 1725
516     beginAtomIndex = 0;
517     beginRigidBodyIndex = 0;
518     beginCutoffGroupIndex = 0;
519    
520 tim 1733 for(int i = 0; i < nproc; i++) {
521 tim 1725 beginAtomIndex += NumAtomsInProc[i];
522     beginRigidBodyIndex += NumRigidBodiesInProc[i];
523 tim 1733 beginCutoffGroupIndex += NumCutoffGroupsInProc[i];
524 tim 1725 }
525 tim 1733
526 tim 1725 #endif
527    
528 tim 1733 for(mol = info->beginMolecule(mi); mol != NULL;
529     mol = info->nextMolecule(mi)) {
530    
531 tim 1725 //local index(index in DataStorge) of atom is important
532 tim 1733 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
533 tim 1725 atom->setGlobalIndex(beginAtomIndex++);
534     }
535    
536 tim 1733 for(rb = mol->beginRigidBody(ri); rb != NULL;
537     rb = mol->nextRigidBody(ri)) {
538 tim 1725 rb->setGlobalIndex(beginRigidBodyIndex++);
539     }
540 tim 1733
541 tim 1725 //local index of cutoff group is trivial, it only depends on the order of travesing
542 tim 1733 for(cg = mol->beginCutoffGroup(ci); cg != NULL;
543     cg = mol->nextCutoffGroup(ci)) {
544     cg->setGlobalIndex(beginCutoffGroupIndex++);
545     }
546     }
547    
548 tim 1735 //fill globalGroupMembership
549 tim 1733 std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
550 tim 1807 for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
551 tim 1725 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
552 tim 1733
553     for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
554     globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
555     }
556    
557     }
558     }
559 tim 1735
560     #ifdef IS_MPI
561 tim 1733 // Since the globalGroupMembership has been zero filled and we've only
562     // poked values into the atoms we know, we can do an Allreduce
563     // to get the full globalGroupMembership array (We think).
564     // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
565     // docs said we could.
566    
567     MPI_Allreduce(&globalGroupMembership[0],
568     info->getGlobalGroupMembershipPointer(),
569     info->getNGlobalAtoms(),
570     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
571 tim 1735 #else
572     std::copy(globalGroupMembership.begin(), globalGroupMembership.end(),
573     info->getGlobalGroupMembershipPointer());
574     #endif
575 tim 1733
576     //fill molMembership
577     std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0);
578    
579     for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
580    
581     for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
582     globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
583     }
584    
585 tim 1735 #ifdef IS_MPI
586 tim 1733 MPI_Allreduce(&globalMolMembership[0],
587     info->getGlobalMolMembershipPointer(),
588     info->getNGlobalAtoms(),
589     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
590 tim 1735 #else
591     std::copy(globalMolMembership.begin(), globalMolMembership.end(),
592     info->getGlobalMolMembershipPointer());
593     #endif
594    
595 tim 1807 }
596 tim 1733 }
597    
598     void SimCreator::loadCoordinates(SimInfo* info) {
599     Globals* globals;
600     globals = info->getGlobals();
601    
602     if (!globals->haveInitialConfig()) {
603     sprintf(painCave.errMsg,
604     "Cannot intialize a simulation without an initial configuration file.\n");
605     painCave.isFatal = 1;;
606     simError();
607     }
608 tim 1725
609 tim 1807 DumpReader reader(info, globals->getInitialConfig());
610     int nframes = reader.getNFrames();
611 tim 1733
612     if (nframes > 0) {
613 tim 1807 reader.readFrame(nframes - 1);
614 tim 1733 } else {
615 tim 1735 //invalid initial coordinate file
616 tim 1733 sprintf(painCave.errMsg, "Initial configuration file %s should at least contain one frame\n",
617     globals->getInitialConfig());
618     painCave.isFatal = 1;
619     simError();
620 tim 1725 }
621 tim 1733
622 tim 1725 }
623    
624 tim 1703 } //end namespace oopse