1 /**
2 * @file Reaction.cpp
3 * @brief Implementations of Reaction and ListOfReactions.
4 * @author Ben Bornstein
5 *
6 * <!--------------------------------------------------------------------------
7 * This file is part of libSBML. Please visit http://sbml.org for more
8 * information about SBML, and the latest version of libSBML.
9 *
10 * Copyright (C) 2020 jointly by the following organizations:
11 * 1. California Institute of Technology, Pasadena, CA, USA
12 * 2. University of Heidelberg, Heidelberg, Germany
13 * 3. University College London, London, UK
14 *
15 * Copyright (C) 2019 jointly by the following organizations:
16 * 1. California Institute of Technology, Pasadena, CA, USA
17 * 2. University of Heidelberg, Heidelberg, Germany
18 *
19 * Copyright (C) 2013-2018 jointly by the following organizations:
20 * 1. California Institute of Technology, Pasadena, CA, USA
21 * 2. EMBL European Bioinformatics Institute (EMBL-EBI), Hinxton, UK
22 * 3. University of Heidelberg, Heidelberg, Germany
23 *
24 * Copyright (C) 2009-2013 jointly by the following organizations:
25 * 1. California Institute of Technology, Pasadena, CA, USA
26 * 2. EMBL European Bioinformatics Institute (EMBL-EBI), Hinxton, UK
27 *
28 * Copyright (C) 2006-2008 by the California Institute of Technology,
29 * Pasadena, CA, USA
30 *
31 * Copyright (C) 2002-2005 jointly by the following organizations:
32 * 1. California Institute of Technology, Pasadena, CA, USA
33 * 2. Japan Science and Technology Agency, Japan
34 *
35 * This library is free software; you can redistribute it and/or modify it
36 * under the terms of the GNU Lesser General Public License as published by
37 * the Free Software Foundation. A copy of the license agreement is provided
38 * in the file named "LICENSE.txt" included with this software distribution
39 * and also available online as http://sbml.org/software/libsbml/license.html
40 * ---------------------------------------------------------------------- -->*/
41
42 #include <sbml/xml/XMLNode.h>
43 #include <sbml/xml/XMLAttributes.h>
44 #include <sbml/xml/XMLInputStream.h>
45 #include <sbml/xml/XMLOutputStream.h>
46
47 #include <sbml/SBO.h>
48 #include <sbml/SBMLVisitor.h>
49 #include <sbml/SpeciesReference.h>
50 #include <sbml/SimpleSpeciesReference.h>
51 #include <sbml/ModifierSpeciesReference.h>
52 #include <sbml/KineticLaw.h>
53 #include <sbml/SBMLDocument.h>
54 #include <sbml/SBMLError.h>
55 #include <sbml/Model.h>
56 #include <sbml/Reaction.h>
57
58 #include <sbml/util/ElementFilter.h>
59
60 /** @cond doxygenIgnored */
61 using namespace std;
62 /** @endcond */
63
64 LIBSBML_CPP_NAMESPACE_BEGIN
65 #ifdef __cplusplus
66
67 /**
68 * Used by getReactant(species), getProduct(species), and
69 * getModifier(species).
70 */
71 static SBase*
GetSpeciesRef(ListOf & items,const string & species)72 GetSpeciesRef (ListOf& items, const string& species)
73 {
74 // TODO: Maybe ListOf should return begin and end iterators to the
75 // underlying container. Then this loop could be rewritten with
76 // a find_if() algorithm.
77
78 unsigned int size = items.size();
79
80 for (unsigned int n = 0; n < size; ++n)
81 {
82 SpeciesReference* sr = static_cast<SpeciesReference*>( items.get(n) );
83 if (sr->getSpecies() == species) return sr;
84 else if (sr->getId() == species) return sr;
85 }
86
87 return NULL;
88 }
89
90
91 /**
92 * Simply calls non-const version above.
93 */
94 static const SBase*
GetSpeciesRef(const ListOf & items,const std::string & species)95 GetSpeciesRef (const ListOf& items, const std::string& species)
96 {
97 return GetSpeciesRef(const_cast<ListOf&>(items), species);
98 }
99
100
101
Reaction(unsigned int level,unsigned int version)102 Reaction::Reaction (unsigned int level, unsigned int version) :
103 SBase ( level, version )
104 , mReactants (level, version)
105 , mProducts (level, version)
106 , mModifiers (level, version)
107 , mKineticLaw( NULL )
108 , mReversible( true )
109 , mFast ( false )
110 , mIsSetFast ( false )
111 , mCompartment ( "" )
112 , mIsSetReversible (false)
113 , mExplicitlySetReversible (false)
114 , mExplicitlySetFast (false)
115 {
116 if (!hasValidLevelVersionNamespaceCombination())
117 throw SBMLConstructorException();
118
119 mReactants.setType( ListOfSpeciesReferences::Reactant );
120 mProducts .setType( ListOfSpeciesReferences::Product );
121 mModifiers.setType( ListOfSpeciesReferences::Modifier );
122 // before level 3 reversible and fast was set by default
123 if (level < 3)
124 {
125 /* this changes existing behaviour as isSetFast already existed in L2 */
126 // mIsSetFast = true;
127 mIsSetReversible = true;
128 }
129
130 connectToChild();
131 }
132
133
Reaction(SBMLNamespaces * sbmlns)134 Reaction::Reaction (SBMLNamespaces * sbmlns) :
135 SBase ( sbmlns )
136 , mReactants (sbmlns)
137 , mProducts (sbmlns)
138 , mModifiers (sbmlns)
139 , mKineticLaw( NULL )
140 , mReversible( true )
141 , mFast ( false )
142 , mIsSetFast ( false )
143 , mCompartment ( "" )
144 , mIsSetReversible (false)
145 , mExplicitlySetReversible (false)
146 , mExplicitlySetFast (false)
147 {
148 if (!hasValidLevelVersionNamespaceCombination())
149 {
150 throw SBMLConstructorException(getElementName(), sbmlns);
151 }
152
153 mReactants.setType( ListOfSpeciesReferences::Reactant );
154 mProducts .setType( ListOfSpeciesReferences::Product );
155 mModifiers.setType( ListOfSpeciesReferences::Modifier );
156
157 // before level 3 reversible and fast was set by default
158 if (sbmlns->getLevel() < 3)
159 {
160 /* this changes existing behaviour as isSetFast already existed in L2 */
161 // mIsSetFast = true;
162 mIsSetReversible = true;
163 }
164 connectToChild();
165 loadPlugins(sbmlns);
166 }
167
168
169 /*
170 * Destroys this Reaction.
171 */
~Reaction()172 Reaction::~Reaction ()
173 {
174 delete mKineticLaw;
175 }
176
177
178 /*
179 * Copy constructor. Creates a copy of this Reaction.
180 */
Reaction(const Reaction & orig)181 Reaction::Reaction (const Reaction& orig)
182 : SBase ( orig )
183 , mReactants ( orig.mReactants )
184 , mProducts ( orig.mProducts )
185 , mModifiers ( orig.mModifiers )
186 , mKineticLaw( NULL )
187 , mReversible( orig.mReversible )
188 , mFast ( orig.mFast )
189 , mIsSetFast ( orig.mIsSetFast )
190 , mCompartment ( orig.mCompartment )
191 , mIsSetReversible ( orig.mIsSetReversible )
192 , mExplicitlySetReversible ( orig.mExplicitlySetReversible )
193 , mExplicitlySetFast ( orig.mExplicitlySetFast )
194 {
195 if (orig.mKineticLaw != NULL)
196 {
197 mKineticLaw = static_cast<KineticLaw*>( orig.mKineticLaw->clone() );
198 }
199 connectToChild();
200 }
201
202
203 /*
204 * Assignment operator.
205 */
operator =(const Reaction & rhs)206 Reaction& Reaction::operator=(const Reaction& rhs)
207 {
208 if(&rhs!=this)
209 {
210 this->SBase::operator =(rhs);
211 mReversible = rhs.mReversible ;
212 mFast = rhs.mFast ;
213 mIsSetFast = rhs.mIsSetFast ;
214 mReactants = rhs.mReactants ;
215 mProducts = rhs.mProducts ;
216 mModifiers = rhs.mModifiers ;
217 mCompartment = rhs.mCompartment;
218 mIsSetReversible = rhs.mIsSetReversible;
219 mExplicitlySetReversible = rhs.mExplicitlySetReversible;
220 mExplicitlySetFast = rhs.mExplicitlySetFast;
221
222 delete mKineticLaw;
223 if (rhs.mKineticLaw != NULL)
224 {
225 mKineticLaw = static_cast<KineticLaw*>( rhs.mKineticLaw->clone() );
226 }
227 else
228 {
229 mKineticLaw = NULL;
230 }
231 }
232
233 connectToChild();
234
235 return *this;
236 }
237
238
239 /** @cond doxygenLibsbmlInternal */
240 bool
accept(SBMLVisitor & v) const241 Reaction::accept (SBMLVisitor& v) const
242 {
243 bool result = v.visit(*this);
244
245 mReactants.accept(v);
246 mProducts .accept(v);
247 mModifiers.accept(v);
248
249 if (mKineticLaw != NULL) mKineticLaw->accept(v);
250
251 v.leave(*this);
252
253 return result;
254 }
255 /** @endcond */
256
257
258 /*
259 * @return a (deep) copy of this Reaction.
260 */
261 Reaction*
clone() const262 Reaction::clone () const
263 {
264 return new Reaction(*this);
265 }
266
267
268 SBase*
getElementBySId(const std::string & id)269 Reaction::getElementBySId(const std::string& id)
270 {
271 if (id.empty()) return NULL;
272 if (mReactants.getId() == id) return &mReactants;
273 if (mProducts.getId() == id) return &mProducts;
274 if (mModifiers.getId() == id) return &mModifiers;
275 if (mKineticLaw && mKineticLaw->getId() == id) return mKineticLaw;
276
277 SBase* obj = mReactants.getElementBySId(id);
278 if (obj != NULL) return obj;
279 obj = mProducts.getElementBySId(id);
280 if (obj != NULL) return obj;
281 obj = mModifiers.getElementBySId(id);
282 if (obj != NULL) return obj;
283 if (mKineticLaw != NULL) {
284 obj = mKineticLaw->getElementBySId(id);
285 if (obj != NULL) return obj;
286 }
287 return getElementFromPluginsBySId(id);
288 }
289
290
291 SBase*
getElementByMetaId(const std::string & metaid)292 Reaction::getElementByMetaId(const std::string& metaid)
293 {
294 if (metaid.empty()) return NULL;
295 if (mReactants.getMetaId() == metaid) return &mReactants;
296 if (mProducts.getMetaId() == metaid) return &mProducts;
297 if (mModifiers.getMetaId() == metaid) return &mModifiers;
298 if (mKineticLaw && mKineticLaw->getMetaId() == metaid) return mKineticLaw;
299
300 SBase* obj = mReactants.getElementByMetaId(metaid);
301 if (obj != NULL) return obj;
302 obj = mProducts.getElementByMetaId(metaid);
303 if (obj != NULL) return obj;
304 obj = mModifiers.getElementByMetaId(metaid);
305 if (obj != NULL) return obj;
306 if (mKineticLaw != NULL) {
307 obj = mKineticLaw->getElementByMetaId(metaid);
308 if (obj != NULL) return obj;
309 }
310 return getElementFromPluginsByMetaId(metaid);
311 }
312
313
314 List*
getAllElements(ElementFilter * filter)315 Reaction::getAllElements(ElementFilter *filter)
316 {
317 List* ret = new List();
318 List* sublist = NULL;
319
320 ADD_FILTERED_POINTER(ret, sublist, mKineticLaw, filter);
321
322 ADD_FILTERED_LIST(ret, sublist, mReactants, filter);
323 ADD_FILTERED_LIST(ret, sublist, mProducts, filter);
324 ADD_FILTERED_LIST(ret, sublist, mModifiers, filter);
325
326 ADD_FILTERED_FROM_PLUGIN(ret, sublist, filter);
327
328 return ret;
329 }
330
331 void
renameSIdRefs(const std::string & oldid,const std::string & newid)332 Reaction::renameSIdRefs(const std::string& oldid, const std::string& newid)
333 {
334 SBase::renameSIdRefs(oldid, newid);
335 if (mCompartment == oldid) {
336 setCompartment(newid);
337 }
338 }
339
340 /*
341 * Initializes the fields of this Reaction to their defaults:
342 *
343 * - reversible = true
344 * - fast = false (L1 only)
345 */
346 void
initDefaults()347 Reaction::initDefaults ()
348 {
349 setReversible(true);
350 // not explicilty set
351 mExplicitlySetReversible = false;
352
353 //
354 // Set fast explicitly and make sure mIsSetFast is false. This preserves
355 // backward compatibility with L1 where fast defaulted to false and such
356 // Reaction.isSetFast() was not available. E.g.:
357 //
358 // Level 1 Level 2
359 // --------------------------- -------------------------------
360 // Reaction r; Reaction r;
361 // r.getFast() == false; r.getFast() == false, but
362 // r.isSetFast() == N/A r.isSetFast() == false
363 //
364 // but for L3 acknowledge that fast has been set
365 mFast = false;
366 mIsSetFast = false;
367 if (getLevel() == 3)
368 {
369 setFast(false);
370 }
371
372 mExplicitlySetFast = false;
373 }
374
375
376 /*
377 * @return the id of this SBML object.
378 */
379 const string&
getId() const380 Reaction::getId () const
381 {
382 return mId;
383 }
384
385
386 /*
387 * @return the name of this SBML object.
388 */
389 const string&
getName() const390 Reaction::getName () const
391 {
392 return (getLevel() == 1) ? mId : mName;
393 }
394
395
396 /*
397 * @return the KineticLaw of this Reaction.
398 */
399 const KineticLaw*
getKineticLaw() const400 Reaction::getKineticLaw () const
401 {
402 return mKineticLaw;
403 }
404
405
406 /*
407 * @return the KineticLaw of this Reaction.
408 */
409 KineticLaw*
getKineticLaw()410 Reaction::getKineticLaw ()
411 {
412 return mKineticLaw;
413 }
414
415
416 /*
417 * @return the reversible status of this Reaction.
418 */
419 bool
getReversible() const420 Reaction::getReversible () const
421 {
422 return mReversible;
423 }
424
425
426 /*
427 * @return the fast status of this Reaction.
428 */
429 bool
getFast() const430 Reaction::getFast () const
431 {
432 return mFast;
433 }
434
435
436 /*
437 * @return the compartment of this SBML object.
438 */
439 const string&
getCompartment() const440 Reaction::getCompartment () const
441 {
442 return mCompartment;
443 }
444
445
446 /*
447 * @return @c true if the id of this SBML object is set, false
448 * otherwise.
449 */
450 bool
isSetId() const451 Reaction::isSetId () const
452 {
453 return (mId.empty() == false);
454 }
455
456
457 /*
458 * @return @c true if the name of this SBML object is set, false
459 * otherwise.
460 */
461 bool
isSetName() const462 Reaction::isSetName () const
463 {
464 return (getLevel() == 1) ? (mId.empty() == false) :
465 (mName.empty() == false);
466 }
467
468
469 /*
470 * @return @c true if the KineticLaw of this Reaction is set, false
471 * otherwise.
472 */
473 bool
isSetKineticLaw() const474 Reaction::isSetKineticLaw () const
475 {
476 return (mKineticLaw != NULL);
477 }
478
479
480 /*
481 * @return @c true if the fast status of this Reaction is set, false
482 * otherwise.
483 *
484 * In L1, fast is optional with a default of false, which means it is
485 * effectively always set. In L2, however, fast is optional with no
486 * default value, so it may or may not be set to a specific value.
487 */
488 bool
isSetFast() const489 Reaction::isSetFast () const
490 {
491 return mIsSetFast;
492 }
493
494
495 /*
496 * @return @c true if the compartment of this SBML object is set, false
497 * otherwise.
498 */
499 bool
isSetCompartment() const500 Reaction::isSetCompartment () const
501 {
502 return (mCompartment.empty() == false);
503 }
504
505
506 /*
507 * @return @c true if the fast status of this Reaction is set, false
508 * otherwise.
509 */
510 bool
isSetReversible() const511 Reaction::isSetReversible () const
512 {
513 return mIsSetReversible;
514 }
515
516
517 /*
518 * Sets the id of this SBML object to a copy of @p sid.
519 */
520 int
setId(const std::string & sid)521 Reaction::setId (const std::string& sid)
522 {
523 /* since the setId function has been used as an
524 * alias for setName we cant require it to only
525 * be used on a L2 model
526 */
527 /* if (getLevel() == 1)
528 {
529 return LIBSBML_UNEXPECTED_ATTRIBUTE;
530 }
531 */
532 if (!(SyntaxChecker::isValidInternalSId(sid)))
533 {
534 return LIBSBML_INVALID_ATTRIBUTE_VALUE;
535 }
536 else
537 {
538 mId = sid;
539 return LIBSBML_OPERATION_SUCCESS;
540 }
541 }
542
543
544 /*
545 * Sets the name of this SBML object to a copy of name.
546 */
547 int
setName(const std::string & name)548 Reaction::setName (const std::string& name)
549 {
550 /* if this is setting an L2 name the type is string
551 * whereas if it is setting an L1 name its type is SId
552 */
553 if (getLevel() == 1)
554 {
555 if (!(SyntaxChecker::isValidInternalSId(name)))
556 {
557 return LIBSBML_INVALID_ATTRIBUTE_VALUE;
558 }
559 else
560 {
561 mId = name;
562 return LIBSBML_OPERATION_SUCCESS;
563 }
564 }
565 else
566 {
567 mName = name;
568 return LIBSBML_OPERATION_SUCCESS;
569 }
570 }
571
572
573 /*
574 * Sets the KineticLaw of this Reaction to a copy of the given KineticLaw.
575 */
576 int
setKineticLaw(const KineticLaw * kl)577 Reaction::setKineticLaw (const KineticLaw* kl)
578 {
579 int returnValue = checkCompatibility(static_cast<const SBase *>(kl));
580
581 if (returnValue == LIBSBML_OPERATION_FAILED && kl == NULL)
582 {
583 delete mKineticLaw;
584 mKineticLaw = NULL;
585 return LIBSBML_OPERATION_SUCCESS;
586 }
587 else if (returnValue != LIBSBML_OPERATION_SUCCESS)
588 {
589 return returnValue;
590 }
591
592 if (mKineticLaw == kl)
593 {
594 return LIBSBML_OPERATION_SUCCESS;
595 }
596 else
597 {
598 delete mKineticLaw;
599 mKineticLaw = static_cast<KineticLaw*>( kl->clone() );
600
601 if (mKineticLaw != NULL) mKineticLaw->connectToParent(this);
602
603 return LIBSBML_OPERATION_SUCCESS;
604 }
605 }
606
607
608 /*
609 * Sets the reversible status of this Reaction to value.
610 */
611 int
setReversible(bool value)612 Reaction::setReversible (bool value)
613 {
614 mReversible = value;
615 mIsSetReversible = true;
616 mExplicitlySetReversible = true;
617 return LIBSBML_OPERATION_SUCCESS;
618 }
619
620
621 /*
622 * Sets the fast status of this Reaction to value.
623 */
624 int
setFast(bool value)625 Reaction::setFast (bool value)
626 {
627 if (getLevel() == 3 && getVersion() > 1)
628 {
629 mFast = false;
630 mIsSetFast = false;
631 mExplicitlySetFast = false;
632 return LIBSBML_UNEXPECTED_ATTRIBUTE;
633 }
634 else
635 {
636 mFast = value;
637 mIsSetFast = true;
638 mExplicitlySetFast = true;
639 return LIBSBML_OPERATION_SUCCESS;
640 }
641 }
642
643
644 /*
645 * Sets the compartment of this SBML object to a copy of @p sid.
646 */
647 int
setCompartment(const std::string & sid)648 Reaction::setCompartment (const std::string& sid)
649 {
650 if (getLevel() < 3)
651 {
652 return LIBSBML_UNEXPECTED_ATTRIBUTE;
653 }
654
655 if (!(SyntaxChecker::isValidInternalSId(sid)))
656 {
657 return LIBSBML_INVALID_ATTRIBUTE_VALUE;
658 }
659 else
660 {
661 mCompartment = sid;
662 return LIBSBML_OPERATION_SUCCESS;
663 }
664 }
665
666
667 /*
668 * Unsets the name of this SBML object.
669 */
670 int
unsetName()671 Reaction::unsetName ()
672 {
673 if (getLevel() == 1)
674 {
675 mId.erase();
676 }
677 else
678 {
679 mName.erase();
680 }
681
682 if (getLevel() == 1 && mId.empty())
683 {
684 return LIBSBML_OPERATION_SUCCESS;
685 }
686 else if (mName.empty())
687 {
688 return LIBSBML_OPERATION_SUCCESS;
689 }
690 else
691 {
692 return LIBSBML_OPERATION_FAILED;
693 }
694 }
695
696
697 /*
698 * Unsets the KineticLaw of this Reaction.
699 */
700 int
unsetKineticLaw()701 Reaction::unsetKineticLaw ()
702 {
703 delete mKineticLaw;
704 mKineticLaw = NULL;
705
706 if (mKineticLaw == NULL)
707 {
708 return LIBSBML_OPERATION_SUCCESS;
709 }
710 else
711 {
712 return LIBSBML_OPERATION_FAILED;
713 }
714 }
715
716
717 /*
718 * Unsets the fast status of this Reaction.
719 *
720 * In L1, fast is optional with a default of false, which means it is
721 * effectively always set. In L2, however, fast is optional with no
722 * default value, so it may or may not be set to a specific value.
723 */
724 int
unsetFast()725 Reaction::unsetFast ()
726 {
727 mIsSetFast = false;
728
729 if (getLevel() == 3 && getVersion() > 1)
730 {
731 return LIBSBML_UNEXPECTED_ATTRIBUTE;
732 }
733
734 if (!mIsSetFast)
735 {
736 return LIBSBML_OPERATION_SUCCESS;
737 }
738 else
739 {
740 return LIBSBML_OPERATION_FAILED;
741 }
742 }
743
744
745 /*
746 * Unsets the compartment of this SBML object.
747 */
748 int
unsetCompartment()749 Reaction::unsetCompartment ()
750 {
751 if (getLevel() < 3)
752 {
753 mCompartment.erase();
754 return LIBSBML_UNEXPECTED_ATTRIBUTE;
755 }
756
757 mCompartment.erase();
758
759 if (mCompartment.empty())
760 {
761 return LIBSBML_OPERATION_SUCCESS;
762 }
763 else
764 {
765 return LIBSBML_OPERATION_FAILED;
766 }
767 }
768
769
770 /*
771 * Unsets the reversible status of this Reaction.
772 */
773 int
unsetReversible()774 Reaction::unsetReversible ()
775 {
776 if (getLevel() < 3)
777 {
778 // reset default
779 mReversible = true;
780 mIsSetReversible = true;
781 mExplicitlySetReversible = false;
782 return LIBSBML_UNEXPECTED_ATTRIBUTE;
783 }
784 else
785 {
786 mIsSetReversible = false;
787 mExplicitlySetReversible = false;
788 return LIBSBML_OPERATION_SUCCESS;
789 }
790 }
791
792
793 /*
794 * Adds a copy of the given reactant (SpeciesReference) to this Reaction.
795 */
796 int
addReactant(const SpeciesReference * sr)797 Reaction::addReactant (const SpeciesReference* sr)
798 {
799 if (sr == NULL)
800 return LIBSBML_OPERATION_FAILED;
801
802 int returnValue = checkCompatibility(static_cast<const SBase *>(sr));
803 if (returnValue != LIBSBML_OPERATION_SUCCESS)
804 {
805 return returnValue;
806 }
807 else if (sr->isSetId()
808 && (getListOfReactants()->get(sr->getId())) != NULL)
809 {
810 // an object with this id already exists
811 return LIBSBML_DUPLICATE_OBJECT_ID;
812 }
813 else
814 {
815 return mReactants.append(sr);
816 }
817 }
818
addReactant(const Species * species,double stoichiometry,const std::string id,bool constant)819 int Reaction::addReactant(
820 const Species *species,
821 double stoichiometry,
822 const std::string id,
823 bool constant)
824 {
825 if (species == NULL)
826 return LIBSBML_INVALID_OBJECT;
827
828 if (!species->isSetId())
829 return LIBSBML_INVALID_ATTRIBUTE_VALUE;
830
831 if (!id.empty() && getListOfReactants()->get(id) != NULL)
832 {
833 // an object with this id already exists
834 return LIBSBML_DUPLICATE_OBJECT_ID;
835 }
836
837 SpeciesReference* ref = createReactant();
838 if (!id.empty())
839 ref->setId(id);
840
841 if (stoichiometry == stoichiometry)
842 ref->setStoichiometry(stoichiometry);
843
844 ref->setSpecies(species->getId());
845 ref->setConstant(constant);
846
847 return LIBSBML_OPERATION_SUCCESS;
848 }
849
850
851 /*
852 * Adds a copy of the given product (SpeciesReference) to this Reaction.
853 */
854 int
addProduct(const SpeciesReference * sr)855 Reaction::addProduct (const SpeciesReference* sr)
856 {
857 if (sr == NULL)
858 return LIBSBML_OPERATION_FAILED;
859
860 int returnValue = checkCompatibility(static_cast<const SBase *>(sr));
861 if (returnValue != LIBSBML_OPERATION_SUCCESS)
862 {
863 return returnValue;
864 }
865 else if (sr->isSetId()
866 && (getListOfProducts()->get(sr->getId())) != NULL)
867 {
868 // an object with this id already exists
869 return LIBSBML_DUPLICATE_OBJECT_ID;
870 }
871 else
872 {
873 return mProducts.append(sr);
874 }
875 }
876
877 int
addProduct(const Species * species,double stoichiometry,const std::string id,bool constant)878 Reaction::addProduct(
879 const Species *species,
880 double stoichiometry,
881 const std::string id,
882 bool constant)
883 {
884 if (species == NULL)
885 return LIBSBML_INVALID_OBJECT;
886
887 if (!species->isSetId())
888 return LIBSBML_INVALID_ATTRIBUTE_VALUE;
889
890 if (!id.empty() && getListOfProducts()->get(id) != NULL)
891 {
892 // an object with this id already exists
893 return LIBSBML_DUPLICATE_OBJECT_ID;
894 }
895
896 SpeciesReference* ref = createProduct();
897 if (!id.empty())
898 ref->setId(id);
899
900 if (stoichiometry == stoichiometry)
901 ref->setStoichiometry(stoichiometry);
902
903 ref->setSpecies(species->getId());
904 ref->setConstant(constant);
905
906 return LIBSBML_OPERATION_SUCCESS;
907 }
908
909
910 /*
911 * Adds a copy of the given modifier (ModifierSpeciesReference) to this
912 * Reaction.
913 */
914 int
addModifier(const ModifierSpeciesReference * msr)915 Reaction::addModifier (const ModifierSpeciesReference* msr)
916 {
917 if (msr == NULL)
918 return LIBSBML_OPERATION_FAILED;
919
920 int returnValue = checkCompatibility(static_cast<const SBase *>(msr));
921 if (returnValue != LIBSBML_OPERATION_SUCCESS)
922 {
923 return returnValue;
924 }
925 else if (msr->isSetId()
926 && (getListOfModifiers()->get(msr->getId())) != NULL)
927 {
928 // an object with this id already exists
929 return LIBSBML_DUPLICATE_OBJECT_ID;
930 }
931 else
932 {
933 return mModifiers.append(msr);
934 }
935 }
936
937 int
addModifier(const Species * species,const std::string id)938 Reaction::addModifier(
939 const Species *species,
940 const std::string id)
941 {
942 if (species == NULL)
943 return LIBSBML_INVALID_OBJECT;
944
945 if (!species->isSetId())
946 return LIBSBML_INVALID_ATTRIBUTE_VALUE;
947
948 if (!id.empty() && getListOfModifiers()->get(id) != NULL)
949 {
950 // an object with this id already exists
951 return LIBSBML_DUPLICATE_OBJECT_ID;
952 }
953
954 ModifierSpeciesReference* ref = createModifier();
955 if (!id.empty())
956 ref->setId(id);
957
958 ref->setSpecies(species->getId());
959
960 return LIBSBML_OPERATION_SUCCESS;
961 }
962
963
964 /*
965 * Creates a new SpeciesReference, adds it to this Reaction's list of
966 * reactants and returns it.
967 */
968 SpeciesReference*
createReactant()969 Reaction::createReactant ()
970 {
971 SpeciesReference* sr = NULL;
972
973 try
974 {
975 sr = new SpeciesReference(getSBMLNamespaces());
976 }
977 catch (...)
978 {
979 /* here we do not create a default object as the level/version must
980 * match the parent object
981 *
982 * so do nothing
983 */
984 return NULL;
985 }
986
987 if (sr != NULL) mReactants.appendAndOwn(sr);
988
989 return sr;
990 }
991
992
993 /*
994 * Creates a new SpeciesReference, adds it to this Reaction's list of
995 * products and returns it.
996 */
997 SpeciesReference*
createProduct()998 Reaction::createProduct ()
999 {
1000 SpeciesReference* sr = NULL;
1001
1002 try
1003 {
1004 sr = new SpeciesReference(getSBMLNamespaces());
1005 }
1006 catch (...)
1007 {
1008 /* here we do not create a default object as the level/version must
1009 * match the parent object
1010 *
1011 * so do nothing
1012 */
1013 return NULL;
1014 }
1015
1016 if (sr) mProducts.appendAndOwn(sr);
1017
1018 return sr;
1019 }
1020
1021
1022 /*
1023 * Creates a new ModifierSpeciesReference, adds it to this Reaction's
1024 * list of modifiers and returns it.
1025 */
1026 ModifierSpeciesReference*
createModifier()1027 Reaction::createModifier ()
1028 {
1029 ModifierSpeciesReference* sr = NULL;
1030
1031 try
1032 {
1033 sr = new ModifierSpeciesReference(getSBMLNamespaces());
1034 }
1035 catch (...)
1036 {
1037 /* here we do not create a default object as the level/version must
1038 * match the parent object
1039 *
1040 * so do nothing
1041 */
1042 return NULL;
1043 }
1044
1045 if (sr != NULL) mModifiers.appendAndOwn(sr);
1046
1047 return sr;
1048 }
1049
1050
1051 /*
1052 * Creates a new KineticLaw for this Reaction and returns it. If this
1053 * Reaction had a previous KineticLaw, it will be destroyed.
1054 */
1055 KineticLaw*
createKineticLaw()1056 Reaction::createKineticLaw ()
1057 {
1058 delete mKineticLaw;
1059 mKineticLaw = NULL;
1060
1061 try
1062 {
1063 mKineticLaw = new KineticLaw(getSBMLNamespaces());
1064 }
1065 catch (...)
1066 {
1067 /* here we do not create a default object as the level/version must
1068 * match the parent object
1069 *
1070 * so do nothing
1071 */
1072 return NULL;
1073 }
1074
1075 if (mKineticLaw != NULL)
1076 {
1077 mKineticLaw->connectToParent(this);
1078 }
1079
1080 return mKineticLaw;
1081 }
1082
1083
1084 /*
1085 * @return the list of Reactants for this Reaction.
1086 */
1087 const ListOfSpeciesReferences*
getListOfReactants() const1088 Reaction::getListOfReactants () const
1089 {
1090 return &mReactants;
1091 }
1092
1093
1094 /*
1095 * @return the list of Reactants for this Reaction.
1096 */
1097 ListOfSpeciesReferences*
getListOfReactants()1098 Reaction::getListOfReactants ()
1099 {
1100 return &mReactants;
1101 }
1102
1103
1104 /*
1105 * @return the list of Products for this Reaction.
1106 */
1107 const ListOfSpeciesReferences*
getListOfProducts() const1108 Reaction::getListOfProducts () const
1109 {
1110 return &mProducts;
1111 }
1112
1113
1114 /*
1115 * @return the list of Products for this Reaction.
1116 */
1117 ListOfSpeciesReferences*
getListOfProducts()1118 Reaction::getListOfProducts ()
1119 {
1120 return &mProducts;
1121 }
1122
1123
1124 /*
1125 * @return the list of Modifiers for this Reaction.
1126 */
1127 const ListOfSpeciesReferences*
getListOfModifiers() const1128 Reaction::getListOfModifiers () const
1129 {
1130 return &mModifiers;
1131 }
1132
1133
1134 /*
1135 * @return the list of Modifiers for this Reaction.
1136 */
1137 ListOfSpeciesReferences*
getListOfModifiers()1138 Reaction::getListOfModifiers ()
1139 {
1140 return &mModifiers;
1141 }
1142
1143
1144 /*
1145 * @return the nth reactant (SpeciesReference) of this Reaction.
1146 */
1147 const SpeciesReference*
getReactant(unsigned int n) const1148 Reaction::getReactant (unsigned int n) const
1149 {
1150 return static_cast<const SpeciesReference*>( mReactants.get(n) );
1151 }
1152
1153
1154 /*
1155 * @return the nth reactant (SpeciesReference) of this Reaction.
1156 */
1157 SpeciesReference*
getReactant(unsigned int n)1158 Reaction::getReactant (unsigned int n)
1159 {
1160 return static_cast<SpeciesReference*>( mReactants.get(n) );
1161 }
1162
1163
1164 /*
1165 * @return the reactant (SpeciesReference) in this Reaction with the given
1166 * species or @c NULL if no such reactant exists.
1167 */
1168 const SpeciesReference*
getReactant(const std::string & species) const1169 Reaction::getReactant (const std::string& species) const
1170 {
1171 return
1172 static_cast<const SpeciesReference*>( GetSpeciesRef(mReactants, species) );
1173 }
1174
1175
1176 /*
1177 * @return the reactant (SpeciesReference) in this Reaction with the given
1178 * species or @c NULL if no such reactant exists.
1179 */
1180 SpeciesReference*
getReactant(const std::string & species)1181 Reaction::getReactant (const std::string& species)
1182 {
1183 return static_cast<SpeciesReference*>( GetSpeciesRef(mReactants, species) );
1184 }
1185
1186
1187 /*
1188 * @return the nth product (SpeciesReference) of this Reaction.
1189 */
1190 const SpeciesReference*
getProduct(unsigned int n) const1191 Reaction::getProduct (unsigned int n) const
1192 {
1193 return static_cast<const SpeciesReference*>( mProducts.get(n) );
1194 }
1195
1196
1197 /*
1198 * @return the nth product (SpeciesReference) of this Reaction.
1199 */
1200 SpeciesReference*
getProduct(unsigned int n)1201 Reaction::getProduct (unsigned int n)
1202 {
1203 return static_cast<SpeciesReference*>( mProducts.get(n) );
1204 }
1205
1206
1207 /*
1208 * @return the product (SpeciesReference) in this Reaction with the given
1209 * species or @c NULL if no such product exists.
1210 */
1211 const SpeciesReference*
getProduct(const std::string & species) const1212 Reaction::getProduct (const std::string& species) const
1213 {
1214 return
1215 static_cast<const SpeciesReference*>( GetSpeciesRef(mProducts, species) );
1216 }
1217
1218
1219 /*
1220 * @return the product (SpeciesReference) in this Reaction with the given
1221 * species or @c NULL if no such product exists.
1222 */
1223 SpeciesReference*
getProduct(const std::string & species)1224 Reaction::getProduct (const std::string& species)
1225 {
1226 return static_cast<SpeciesReference*>( GetSpeciesRef(mProducts, species) );
1227 }
1228
1229
1230 /*
1231 * @return the nth modifier (ModifierSpeciesReference) of this Reaction.
1232 */
1233 const ModifierSpeciesReference*
getModifier(unsigned int n) const1234 Reaction::getModifier (unsigned int n) const
1235 {
1236 return static_cast<const ModifierSpeciesReference*>( mModifiers.get(n) );
1237 }
1238
1239
1240 /*
1241 * @return the nth modifier (ModifierSpeciesReference) of this Reaction.
1242 */
1243 ModifierSpeciesReference*
getModifier(unsigned int n)1244 Reaction::getModifier (unsigned int n)
1245 {
1246 return static_cast<ModifierSpeciesReference*>( mModifiers.get(n) );
1247 }
1248
1249
1250 /*
1251 * @return the modifier (ModifierSpeciesReference) in this Reaction with
1252 * the given species or @c NULL if no such modifier exists.
1253 */
1254 const ModifierSpeciesReference*
getModifier(const std::string & species) const1255 Reaction::getModifier (const std::string& species) const
1256 {
1257 return static_cast<const ModifierSpeciesReference*>
1258 (
1259 GetSpeciesRef(mModifiers, species)
1260 );
1261 }
1262
1263
1264 /*
1265 * @return the modifier (ModifierSpeciesReference) in this Reaction with
1266 * the given species or @c NULL if no such modifier exists.
1267 */
1268 ModifierSpeciesReference*
getModifier(const std::string & species)1269 Reaction::getModifier (const std::string& species)
1270 {
1271 return static_cast<ModifierSpeciesReference*>
1272 (
1273 GetSpeciesRef(mModifiers, species)
1274 );
1275 }
1276
1277
1278 /*
1279 * @return the number of reactants (SpeciesReferences) in this Reaction.
1280 */
1281 unsigned int
getNumReactants() const1282 Reaction::getNumReactants () const
1283 {
1284 return mReactants.size();
1285 }
1286
1287
1288 /*
1289 * @return the number of products (SpeciesReferences) in this Reaction.
1290 */
1291 unsigned int
getNumProducts() const1292 Reaction::getNumProducts () const
1293 {
1294 return mProducts.size();
1295 }
1296
1297
1298 /*
1299 * @return the number of modifiers (ModifierSpeciesReferences) in this
1300 * Reaction.
1301 */
1302 unsigned int
getNumModifiers() const1303 Reaction::getNumModifiers () const
1304 {
1305 return mModifiers.size();
1306 }
1307
1308
1309 /**
1310 * Removes the nth reactant species (SpeciesReference object) in the list of
1311 * reactants in this Reaction and returns a pointer to it.
1312 */
1313 SpeciesReference*
removeReactant(unsigned int n)1314 Reaction::removeReactant (unsigned int n)
1315 {
1316 return static_cast<SpeciesReference*>(mReactants.remove(n));
1317 }
1318
1319
1320 /**
1321 * Removes the reactant species (SpeciesReference object) having the given
1322 * "species" attribute in this Reaction and returns a pointer to it.
1323 */
1324 SpeciesReference*
removeReactant(const std::string & species)1325 Reaction::removeReactant (const std::string& species)
1326 {
1327 unsigned int size = mReactants.size();
1328
1329 for (unsigned int n = 0; n < size; ++n)
1330 {
1331 SpeciesReference* sr = static_cast<SpeciesReference*>( mReactants.get(n) );
1332 if (sr->getSpecies() == species)
1333 return static_cast<SpeciesReference*>(mReactants.remove(n));
1334 }
1335 return NULL;
1336 }
1337
1338
1339 /**
1340 * Removes the nth product species (SpeciesReference object) in the list of
1341 * products in this Reaction and returns a pointer to it.
1342 */
1343 SpeciesReference*
removeProduct(unsigned int n)1344 Reaction::removeProduct (unsigned int n)
1345 {
1346 return static_cast<SpeciesReference*>(mProducts.remove(n));
1347 }
1348
1349
1350 /**
1351 * Removes the product species (SpeciesReference object) having the given
1352 * "species" attribute in this Reaction and returns a pointer to it.
1353 */
1354 SpeciesReference*
removeProduct(const std::string & species)1355 Reaction::removeProduct (const std::string& species)
1356 {
1357 unsigned int size = mProducts.size();
1358
1359 for (unsigned int n = 0; n < size; ++n)
1360 {
1361 SpeciesReference* sr = static_cast<SpeciesReference*>( mProducts.get(n) );
1362 if (sr->getSpecies() == species)
1363 return static_cast<SpeciesReference*>(mProducts.remove(n));
1364 }
1365 return NULL;
1366 }
1367
1368
1369 /**
1370 * Removes the nth modifier species (ModifierSpeciesReference object) in
1371 * the list of modifiers in this Reaction and returns a pointer to it.
1372 */
1373 ModifierSpeciesReference*
removeModifier(unsigned int n)1374 Reaction::removeModifier (unsigned int n)
1375 {
1376 return static_cast<ModifierSpeciesReference*>(mModifiers.remove(n));
1377 }
1378
1379
1380 /**
1381 * Removes the modifier species (ModifierSpeciesReference object) having
1382 * the given "species" attribute in this Reaction and returns a pointer to it.
1383 */
1384 ModifierSpeciesReference*
removeModifier(const std::string & species)1385 Reaction::removeModifier (const std::string& species)
1386 {
1387 unsigned int size = mModifiers.size();
1388
1389 for (unsigned int n = 0; n < size; ++n)
1390 {
1391 SpeciesReference* sr = static_cast<SpeciesReference*>( mModifiers.get(n) );
1392 if (sr->getSpecies() == species)
1393 return static_cast<ModifierSpeciesReference*>(mModifiers.remove(n));
1394 }
1395 return NULL;
1396 }
1397
1398
1399 /** @cond doxygenLibsbmlInternal */
1400 /*
1401 * Sets the parent SBMLDocument of this SBML object.
1402 */
1403 void
setSBMLDocument(SBMLDocument * d)1404 Reaction::setSBMLDocument (SBMLDocument* d)
1405 {
1406 SBase::setSBMLDocument(d);
1407
1408 mReactants.setSBMLDocument(d);
1409 mProducts .setSBMLDocument(d);
1410 mModifiers.setSBMLDocument(d);
1411
1412 if (mKineticLaw != NULL) mKineticLaw->setSBMLDocument(d);
1413 }
1414
1415
1416 /*
1417 * Sets this SBML object to child SBML objects (if any).
1418 * (Creates a child-parent relationship by the parent)
1419 */
1420 void
connectToChild()1421 Reaction::connectToChild()
1422 {
1423 SBase::connectToChild();
1424 mReactants.connectToParent(this);
1425 mProducts .connectToParent(this);
1426 mModifiers.connectToParent(this);
1427
1428 if (mKineticLaw) mKineticLaw->connectToParent(this);
1429 }
1430
1431
1432 /**
1433 * Enables/Disables the given package with this element and child
1434 * elements (if any).
1435 * (This is an internal implementation for enablePackage function)
1436 */
1437 void
enablePackageInternal(const std::string & pkgURI,const std::string & pkgPrefix,bool flag)1438 Reaction::enablePackageInternal(const std::string& pkgURI,
1439 const std::string& pkgPrefix, bool flag)
1440 {
1441 SBase::enablePackageInternal(pkgURI,pkgPrefix,flag);
1442
1443 mReactants.enablePackageInternal(pkgURI,pkgPrefix,flag);
1444 mProducts.enablePackageInternal(pkgURI,pkgPrefix,flag);
1445 mModifiers.enablePackageInternal(pkgURI,pkgPrefix,flag);
1446
1447 if (mKineticLaw) mKineticLaw->enablePackageInternal(pkgURI,pkgPrefix,flag);
1448 }
1449
1450
1451 void
updateSBMLNamespace(const std::string & pkg,unsigned int level,unsigned int version)1452 Reaction::updateSBMLNamespace(const std::string& pkg, unsigned int level,
1453 unsigned int version)
1454 {
1455 SBase::updateSBMLNamespace(pkg, level, version);
1456
1457 mReactants.updateSBMLNamespace(pkg, level, version);
1458 mProducts.updateSBMLNamespace(pkg, level, version);
1459 mModifiers.updateSBMLNamespace(pkg, level, version);
1460
1461 if (mKineticLaw) mKineticLaw->updateSBMLNamespace(pkg, level, version);
1462 }
1463 /** @endcond */
1464
1465
1466 /*
1467 * @return the typecode (int) of this SBML object or SBML_UNKNOWN
1468 * (default).
1469 *
1470 * @see getElementName()
1471 */
1472 int
getTypeCode() const1473 Reaction::getTypeCode () const
1474 {
1475 return SBML_REACTION;
1476 }
1477
1478
1479 /*
1480 * @return the name of this element ie "reaction".
1481 */
1482 const string&
getElementName() const1483 Reaction::getElementName () const
1484 {
1485 static const string name = "reaction";
1486 return name;
1487 }
1488
1489
1490 bool
hasRequiredAttributes() const1491 Reaction::hasRequiredAttributes() const
1492 {
1493 bool allPresent = true;
1494
1495 /* required attributes for reaction:
1496 * @li id (name in L1)
1497 * @li fast (in L3V1 only)
1498 * @li reversible (in L3 only)
1499 */
1500
1501 if (!isSetId())
1502 allPresent = false;
1503
1504 if (getLevel() > 2 && !isSetReversible())
1505 allPresent = false;
1506
1507 if (getLevel() == 3 && getVersion() == 1 && !isSetFast())
1508 allPresent = false;
1509
1510 return allPresent;
1511 }
1512
1513
1514 /** @cond doxygenLibsbmlInternal */
1515 /*
1516 * @return the SBML object corresponding to next XMLToken in the
1517 * XMLInputStream or @c NULL if the token was not recognized.
1518 */
1519 SBase*
createObject(XMLInputStream & stream)1520 Reaction::createObject (XMLInputStream& stream)
1521 {
1522 const string& name = stream.peek().getName();
1523 SBase* object = NULL;
1524
1525
1526 if (name == "listOfReactants")
1527 {
1528 if (mReactants.size() != 0)
1529 {
1530 if (getLevel() < 3)
1531 logError(NotSchemaConformant);
1532 else
1533 logError(OneSubElementPerReaction, getLevel(), getVersion());
1534 }
1535 mReactants.setExplicitlyListed();
1536 object = &mReactants;
1537 }
1538 else if (name == "listOfProducts")
1539 {
1540 if (mProducts.size() != 0)
1541 {
1542 if (getLevel() < 3)
1543 logError(NotSchemaConformant);
1544 else
1545 logError(OneSubElementPerReaction, getLevel(), getVersion());
1546 }
1547 mProducts.setExplicitlyListed();
1548 object = &mProducts;
1549 }
1550 else if (name == "listOfModifiers")
1551 {
1552 if (getLevel() == 1)
1553 {
1554 return NULL;
1555 }
1556
1557 if (mModifiers.size() != 0)
1558 {
1559 if (getLevel() < 3)
1560 logError(NotSchemaConformant);
1561 else
1562 logError(OneSubElementPerReaction, getLevel(), getVersion());
1563 }
1564 mModifiers.setExplicitlyListed();
1565 object = &mModifiers;
1566 }
1567 else if (name == "kineticLaw")
1568 {
1569 if (mKineticLaw != NULL)
1570 {
1571 if (getLevel() < 3)
1572 logError(NotSchemaConformant);
1573 else
1574 logError(OneSubElementPerReaction, getLevel(), getVersion());
1575 }
1576 delete mKineticLaw;
1577
1578 try
1579 {
1580 mKineticLaw = new KineticLaw(getSBMLNamespaces());
1581 }
1582 catch (SBMLConstructorException*)
1583 {
1584 mKineticLaw = new KineticLaw(SBMLDocument::getDefaultLevel(),
1585 SBMLDocument::getDefaultVersion());
1586 }
1587 catch ( ... )
1588 {
1589 mKineticLaw = new KineticLaw(SBMLDocument::getDefaultLevel(),
1590 SBMLDocument::getDefaultVersion());
1591 }
1592
1593 object = mKineticLaw;
1594 }
1595
1596 return object;
1597 }
1598 /** @endcond */
1599
1600
1601
1602
1603
1604
1605
1606 /** @cond doxygenLibsbmlInternal */
1607
1608 /*
1609 * Returns the value of the "attributeName" attribute of this Reaction.
1610 */
1611 int
getAttribute(const std::string & attributeName,bool & value) const1612 Reaction::getAttribute(const std::string& attributeName, bool& value) const
1613 {
1614 int return_value = SBase::getAttribute(attributeName, value);
1615
1616 if (return_value == LIBSBML_OPERATION_SUCCESS)
1617 {
1618 return return_value;
1619 }
1620
1621 if (attributeName == "fast")
1622 {
1623 value = getFast();
1624 return_value = LIBSBML_OPERATION_SUCCESS;
1625 }
1626 else if (attributeName == "reversible")
1627 {
1628 value = getReversible();
1629 return_value = LIBSBML_OPERATION_SUCCESS;
1630 }
1631
1632 return return_value;
1633 }
1634
1635 /** @endcond */
1636
1637
1638
1639 /** @cond doxygenLibsbmlInternal */
1640
1641 /*
1642 * Returns the value of the "attributeName" attribute of this Reaction.
1643 */
1644 int
getAttribute(const std::string & attributeName,int & value) const1645 Reaction::getAttribute(const std::string& attributeName, int& value) const
1646 {
1647 int return_value = SBase::getAttribute(attributeName, value);
1648
1649 return return_value;
1650 }
1651
1652 /** @endcond */
1653
1654
1655
1656 /** @cond doxygenLibsbmlInternal */
1657
1658 /*
1659 * Returns the value of the "attributeName" attribute of this Reaction.
1660 */
1661 int
getAttribute(const std::string & attributeName,double & value) const1662 Reaction::getAttribute(const std::string& attributeName, double& value) const
1663 {
1664 int return_value = SBase::getAttribute(attributeName, value);
1665
1666 return return_value;
1667 }
1668
1669 /** @endcond */
1670
1671
1672
1673 /** @cond doxygenLibsbmlInternal */
1674
1675 /*
1676 * Returns the value of the "attributeName" attribute of this Reaction.
1677 */
1678 int
getAttribute(const std::string & attributeName,unsigned int & value) const1679 Reaction::getAttribute(const std::string& attributeName,
1680 unsigned int& value) const
1681 {
1682 int return_value = SBase::getAttribute(attributeName, value);
1683
1684 return return_value;
1685 }
1686
1687 /** @endcond */
1688
1689
1690
1691 /** @cond doxygenLibsbmlInternal */
1692
1693 /*
1694 * Returns the value of the "attributeName" attribute of this Reaction.
1695 */
1696 int
getAttribute(const std::string & attributeName,std::string & value) const1697 Reaction::getAttribute(const std::string& attributeName,
1698 std::string& value) const
1699 {
1700 int return_value = SBase::getAttribute(attributeName, value);
1701
1702 if (return_value == LIBSBML_OPERATION_SUCCESS)
1703 {
1704 return return_value;
1705 }
1706
1707 if (attributeName == "compartment")
1708 {
1709 value = getCompartment();
1710 return_value = LIBSBML_OPERATION_SUCCESS;
1711 }
1712
1713 return return_value;
1714 }
1715
1716 /** @endcond */
1717
1718
1719
1720 /** @cond doxygenLibsbmlInternal */
1721
1722 /*
1723 * Returns the value of the "attributeName" attribute of this Reaction.
1724 */
1725 //int
1726 //Reaction::getAttribute(const std::string& attributeName,
1727 // const char* value) const
1728 //{
1729 // int return_value = SBase::getAttribute(attributeName, value);
1730 //
1731 // if (return_value == LIBSBML_OPERATION_SUCCESS)
1732 // {
1733 // return return_value;
1734 // }
1735 //
1736 // if (attributeName == "compartment")
1737 // {
1738 // value = getCompartment().c_str();
1739 // return_value = LIBSBML_OPERATION_SUCCESS;
1740 // }
1741 //
1742 // return return_value;
1743 //}
1744
1745 /** @endcond */
1746
1747
1748
1749 /** @cond doxygenLibsbmlInternal */
1750
1751 /*
1752 * Predicate returning @c true if this Reaction's attribute "attributeName" is
1753 * set.
1754 */
1755 bool
isSetAttribute(const std::string & attributeName) const1756 Reaction::isSetAttribute(const std::string& attributeName) const
1757 {
1758 bool value = SBase::isSetAttribute(attributeName);
1759
1760 if (attributeName == "fast")
1761 {
1762 value = isSetFast();
1763 }
1764 else if (attributeName == "reversible")
1765 {
1766 value = isSetReversible();
1767 }
1768 else if (attributeName == "compartment")
1769 {
1770 value = isSetCompartment();
1771 }
1772
1773 return value;
1774 }
1775
1776 /** @endcond */
1777
1778
1779
1780 /** @cond doxygenLibsbmlInternal */
1781
1782 /*
1783 * Sets the value of the "attributeName" attribute of this Reaction.
1784 */
1785 int
setAttribute(const std::string & attributeName,bool value)1786 Reaction::setAttribute(const std::string& attributeName, bool value)
1787 {
1788 int return_value = SBase::setAttribute(attributeName, value);
1789
1790 if (attributeName == "fast")
1791 {
1792 return_value = setFast(value);
1793 }
1794 else if (attributeName == "reversible")
1795 {
1796 return_value = setReversible(value);
1797 }
1798
1799 return return_value;
1800 }
1801
1802 /** @endcond */
1803
1804
1805
1806 /** @cond doxygenLibsbmlInternal */
1807
1808 /*
1809 * Sets the value of the "attributeName" attribute of this Reaction.
1810 */
1811 int
setAttribute(const std::string & attributeName,int value)1812 Reaction::setAttribute(const std::string& attributeName, int value)
1813 {
1814 int return_value = SBase::setAttribute(attributeName, value);
1815
1816 return return_value;
1817 }
1818
1819 /** @endcond */
1820
1821
1822
1823 /** @cond doxygenLibsbmlInternal */
1824
1825 /*
1826 * Sets the value of the "attributeName" attribute of this Reaction.
1827 */
1828 int
setAttribute(const std::string & attributeName,double value)1829 Reaction::setAttribute(const std::string& attributeName, double value)
1830 {
1831 int return_value = SBase::setAttribute(attributeName, value);
1832
1833 return return_value;
1834 }
1835
1836 /** @endcond */
1837
1838
1839
1840 /** @cond doxygenLibsbmlInternal */
1841
1842 /*
1843 * Sets the value of the "attributeName" attribute of this Reaction.
1844 */
1845 int
setAttribute(const std::string & attributeName,unsigned int value)1846 Reaction::setAttribute(const std::string& attributeName, unsigned int value)
1847 {
1848 int return_value = SBase::setAttribute(attributeName, value);
1849
1850 return return_value;
1851 }
1852
1853 /** @endcond */
1854
1855
1856
1857 /** @cond doxygenLibsbmlInternal */
1858
1859 /*
1860 * Sets the value of the "attributeName" attribute of this Reaction.
1861 */
1862 int
setAttribute(const std::string & attributeName,const std::string & value)1863 Reaction::setAttribute(const std::string& attributeName,
1864 const std::string& value)
1865 {
1866 int return_value = SBase::setAttribute(attributeName, value);
1867
1868 if (attributeName == "compartment")
1869 {
1870 return_value = setCompartment(value);
1871 }
1872
1873 return return_value;
1874 }
1875
1876 /** @endcond */
1877
1878
1879
1880 /** @cond doxygenLibsbmlInternal */
1881
1882 /*
1883 * Sets the value of the "attributeName" attribute of this Reaction.
1884 */
1885 //int
1886 //Reaction::setAttribute(const std::string& attributeName, const char* value)
1887 //{
1888 // int return_value = SBase::setAttribute(attributeName, value);
1889 //
1890 // if (attributeName == "compartment")
1891 // {
1892 // return_value = setCompartment(value);
1893 // }
1894 //
1895 // return return_value;
1896 //}
1897
1898 /** @endcond */
1899
1900
1901
1902 /** @cond doxygenLibsbmlInternal */
1903
1904 /*
1905 * Unsets the value of the "attributeName" attribute of this Reaction.
1906 */
1907 int
unsetAttribute(const std::string & attributeName)1908 Reaction::unsetAttribute(const std::string& attributeName)
1909 {
1910 int value = SBase::unsetAttribute(attributeName);
1911
1912 if (attributeName == "fast")
1913 {
1914 value = unsetFast();
1915 }
1916 else if (attributeName == "reversible")
1917 {
1918 value = unsetReversible();
1919 }
1920 else if (attributeName == "compartment")
1921 {
1922 value = unsetCompartment();
1923 }
1924
1925 return value;
1926 }
1927
1928 /** @endcond */
1929
1930
1931
1932 /** @cond doxygenLibsbmlInternal */
1933
1934 /*
1935 * Creates and returns an new "elementName" object in this Reaction.
1936 */
1937 SBase*
createChildObject(const std::string & elementName)1938 Reaction::createChildObject(const std::string& elementName)
1939 {
1940 SBase* obj = NULL;
1941
1942 if (elementName == "kineticLaw")
1943 {
1944 return createKineticLaw();
1945 }
1946 else if (elementName == "product")
1947 {
1948 return createProduct();
1949 }
1950 else if (elementName == "reactant")
1951 {
1952 return createReactant();
1953 }
1954 else if (elementName == "modifier")
1955 {
1956 return createModifier();
1957 }
1958
1959 return obj;
1960 }
1961
1962 /** @endcond */
1963
1964 /** @cond doxygenLibsbmlInternal */
1965
1966 /*
1967 * Adds an new "elementName" object in this Reaction.
1968 */
1969 int
addChildObject(const std::string & elementName,const SBase * element)1970 Reaction::addChildObject(const std::string& elementName, const SBase* element)
1971 {
1972 if (elementName == "kineticLaw" && element->getTypeCode() == SBML_KINETIC_LAW)
1973 {
1974 return setKineticLaw((const KineticLaw*)(element));
1975 }
1976 else if (elementName == "reactant" && element->getTypeCode() == SBML_SPECIES_REFERENCE)
1977 {
1978 return addReactant((const SpeciesReference*)(element));
1979 }
1980 else if (elementName == "product" && element->getTypeCode() == SBML_SPECIES_REFERENCE)
1981 {
1982 return addProduct((const SpeciesReference*)(element));
1983 }
1984 else if (elementName == "modifier" && element->getTypeCode() == SBML_MODIFIER_SPECIES_REFERENCE)
1985 {
1986 return addModifier((const ModifierSpeciesReference*)(element));
1987 }
1988
1989 return LIBSBML_OPERATION_FAILED;
1990 }
1991
1992 /** @endcond */
1993
1994
1995 /** @cond doxygenLibsbmlInternal */
1996
1997 /*
1998 * Adds an new "elementName" object in this Reaction.
1999 */
2000 SBase*
removeChildObject(const std::string & elementName,const std::string & id)2001 Reaction::removeChildObject(const std::string& elementName, const std::string& id)
2002 {
2003
2004 if (elementName == "kineticLaw")
2005 {
2006 unsetKineticLaw(); // already deletes the kineticLaw, so nothing to do after
2007 }
2008 else if (elementName == "reactant")
2009 {
2010 return removeReactant(id);
2011 }
2012 else if (elementName == "product")
2013 {
2014 return removeProduct(id);
2015 }
2016 else if (elementName == "modifier")
2017 {
2018 return removeModifier(id);
2019 }
2020
2021 return NULL;
2022 }
2023
2024 /** @endcond */
2025
2026
2027
2028 /** @cond doxygenLibsbmlInternal */
2029
2030 /*
2031 * Returns the number of "elementName" in this Reaction.
2032 */
2033 unsigned int
getNumObjects(const std::string & elementName)2034 Reaction::getNumObjects(const std::string& elementName)
2035 {
2036 unsigned int n = 0;
2037
2038 if (elementName == "kineticLaw")
2039 {
2040 if (isSetKineticLaw())
2041 {
2042 return 1;
2043 }
2044 }
2045 else if (elementName == "reactant")
2046 {
2047 return getNumReactants();
2048 }
2049 else if (elementName == "product")
2050 {
2051 return getNumProducts();
2052 }
2053 else if (elementName == "modifier")
2054 {
2055 return getNumModifiers();
2056 }
2057
2058 return n;
2059 }
2060
2061 /** @endcond */
2062
2063
2064
2065 /** @cond doxygenLibsbmlInternal */
2066
2067 /*
2068 * Returns the nth object of "objectName" in this Reaction.
2069 */
2070 SBase*
getObject(const std::string & elementName,unsigned int index)2071 Reaction::getObject(const std::string& elementName, unsigned int index)
2072 {
2073 SBase* obj = NULL;
2074
2075 if (elementName == "kineticLaw")
2076 {
2077 return getKineticLaw();
2078 }
2079 else if (elementName == "reactant")
2080 {
2081 return getReactant(index);
2082 }
2083 else if (elementName == "product")
2084 {
2085 return getProduct(index);
2086 }
2087 else if (elementName == "modifier")
2088 {
2089 return getModifier(index);
2090 }
2091
2092 return obj;
2093 }
2094
2095 /** @endcond */
2096
2097
2098 /** @cond doxygenLibsbmlInternal */
2099 /**
2100 * Subclasses should override this method to get the list of
2101 * expected attributes.
2102 * This function is invoked from corresponding readAttributes()
2103 * function.
2104 */
2105 void
addExpectedAttributes(ExpectedAttributes & attributes)2106 Reaction::addExpectedAttributes(ExpectedAttributes& attributes)
2107 {
2108 SBase::addExpectedAttributes(attributes);
2109
2110 const unsigned int level = getLevel ();
2111 const unsigned int version = getVersion();
2112
2113 switch (level)
2114 {
2115 case 1:
2116 attributes.add("name");
2117 attributes.add("reversible");
2118 attributes.add("fast");
2119 break;
2120 case 2:
2121 attributes.add("name");
2122 attributes.add("reversible");
2123 attributes.add("fast");
2124 attributes.add("id");
2125 if (version == 2)
2126 {
2127 attributes.add("sboTerm");
2128 }
2129
2130 break;
2131 case 3:
2132 attributes.add("reversible");
2133 attributes.add("compartment");
2134 if (version == 1)
2135 {
2136 attributes.add("name");
2137 attributes.add("id");
2138 attributes.add("fast");
2139 }
2140 break;
2141 default:
2142 // assumes l3v2
2143 attributes.add("reversible");
2144 attributes.add("compartment");
2145 break;
2146 }
2147
2148 }
2149 /** @endcond */
2150
2151
2152 /** @cond doxygenLibsbmlInternal */
2153 /*
2154 * Subclasses should override this method to read values from the given
2155 * XMLAttributes set into their specific fields. Be sure to call your
2156 * parent's implementation of this method as well.
2157 */
2158 void
readAttributes(const XMLAttributes & attributes,const ExpectedAttributes & expectedAttributes)2159 Reaction::readAttributes (const XMLAttributes& attributes,
2160 const ExpectedAttributes& expectedAttributes)
2161 {
2162 const unsigned int level = getLevel ();
2163
2164 SBase::readAttributes(attributes, expectedAttributes);
2165
2166 switch (level)
2167 {
2168 case 1:
2169 readL1Attributes(attributes);
2170 break;
2171 case 2:
2172 readL2Attributes(attributes);
2173 break;
2174 case 3:
2175 default:
2176 readL3Attributes(attributes);
2177 break;
2178 }
2179 }
2180 /** @endcond */
2181
2182
2183 /** @cond doxygenLibsbmlInternal */
2184 /*
2185 * Subclasses should override this method to read values from the given
2186 * XMLAttributes set into their specific fields. Be sure to call your
2187 * parent's implementation of this method as well.
2188 */
2189 void
readL1Attributes(const XMLAttributes & attributes)2190 Reaction::readL1Attributes (const XMLAttributes& attributes)
2191 {
2192 const unsigned int level = getLevel ();
2193 const unsigned int version = getVersion();
2194
2195 //
2196 // name: SName { use="required" } (L1v1, L1v2)
2197 //
2198 bool assigned = attributes.readInto("name", mId, getErrorLog(), true, getLine(), getColumn());
2199 if (assigned && mId.size() == 0)
2200 {
2201 logEmptyString("name", level, version, "<reaction>");
2202 }
2203 if (!SyntaxChecker::isValidInternalSId(mId))
2204 logError(InvalidIdSyntax, level, version, "The id '" + mId + "' does not conform to the syntax.");
2205
2206 //
2207 // reversible: boolean { use="optional" default="true" }
2208 // (L1v1, L1v2, L2v1->)
2209 //
2210 mExplicitlySetReversible = attributes.readInto("reversible", mReversible, getErrorLog(), false, getLine(), getColumn());
2211
2212 //
2213 // fast: boolean { use="optional" default="false" } (L1v1, L1v2)
2214 // fast: boolean { use="optional" } (L2v1 ->)
2215 //
2216 mIsSetFast = attributes.readInto("fast", mFast, getErrorLog(), false, getLine(), getColumn());
2217 mExplicitlySetFast = mIsSetFast;
2218 }
2219 /** @endcond */
2220
2221
2222 /** @cond doxygenLibsbmlInternal */
2223 /*
2224 * Subclasses should override this method to read values from the given
2225 * XMLAttributes set into their specific fields. Be sure to call your
2226 * parent's implementation of this method as well.
2227 */
2228 void
readL2Attributes(const XMLAttributes & attributes)2229 Reaction::readL2Attributes (const XMLAttributes& attributes)
2230 {
2231 const unsigned int level = getLevel ();
2232 const unsigned int version = getVersion();
2233
2234 //
2235 // id: SId { use="required" } (L2v1 ->)
2236 //
2237 bool assigned = attributes.readInto("id", mId, getErrorLog(), true, getLine(), getColumn());
2238 if (assigned && mId.size() == 0)
2239 {
2240 logEmptyString("id", level, version, "<reaction>");
2241 }
2242 if (!SyntaxChecker::isValidInternalSId(mId))
2243 logError(InvalidIdSyntax, level, version, "The id '" + mId + "' does not conform to the syntax.");
2244
2245 //
2246 // reversible: boolean { use="optional" default="true" }
2247 // (L1v1, L1v2, L2v1->)
2248 // reversible: boolean { use="required"} (L3v1->)
2249 //
2250 mExplicitlySetReversible = attributes.readInto("reversible", mReversible, getErrorLog(), false, getLine(), getColumn());
2251
2252 //
2253 // fast: boolean { use="optional" } (L2v1 ->)
2254 //
2255 mIsSetFast = attributes.readInto("fast", mFast, getErrorLog(), false, getLine(), getColumn());
2256 mExplicitlySetFast = mIsSetFast;
2257
2258 //
2259 // name: string { use="optional" } (L2v1 ->)
2260 //
2261 attributes.readInto("name", mName, getErrorLog(), false, getLine(), getColumn());
2262
2263 //
2264 // sboTerm: SBOTerm { use="optional" } (L2v2 ->)
2265 //
2266 if (version == 2)
2267 mSBOTerm = SBO::readTerm(attributes, this->getErrorLog(), level, version,
2268 getLine(), getColumn());
2269
2270 }
2271 /** @endcond */
2272
2273
2274 /** @cond doxygenLibsbmlInternal */
2275 /*
2276 * Subclasses should override this method to read values from the given
2277 * XMLAttributes set into their specific fields. Be sure to call your
2278 * parent's implementation of this method as well.
2279 */
2280 void
readL3Attributes(const XMLAttributes & attributes)2281 Reaction::readL3Attributes (const XMLAttributes& attributes)
2282 {
2283 const unsigned int level = getLevel ();
2284 const unsigned int version = getVersion();
2285
2286 //
2287 // id: SId { use="required" } (L2v1 ->)
2288 //
2289 bool assigned;
2290 // for l3v2 sbase will read this as generically optional
2291 // we want to log errors relating to the specific object
2292 if (version == 1)
2293 {
2294 assigned = attributes.readInto("id", mId, getErrorLog(), false, getLine(), getColumn());
2295 if (!assigned)
2296 {
2297 logError(AllowedAttributesOnReaction, level, version,
2298 "The required attribute 'id' is missing.");
2299 }
2300 if (assigned && mId.size() == 0)
2301 {
2302 logEmptyString("id", level, version, "<reaction>");
2303 }
2304 if (!SyntaxChecker::isValidInternalSId(mId))
2305 logError(InvalidIdSyntax, level, version, "The id '" + mId + "' does not conform to the syntax.");
2306 }
2307 else
2308 {
2309 // need to check that id was present
2310 // it has already been read and checked for syntax/emptyness
2311 if (attributes.hasAttribute("id") == false)
2312 {
2313 logError(AllowedAttributesOnReaction, level, version,
2314 "The required attribute 'id' is missing.");
2315 }
2316 }
2317
2318 string elplusid = "<reaction>";
2319 if (!mId.empty()) {
2320 elplusid += " with the id '" + mId + "'";
2321 }
2322 //
2323 // reversible: boolean { use="required"} (L3v1->)
2324 //
2325 mIsSetReversible = attributes.readInto("reversible",
2326 mReversible, getErrorLog(), false, getLine(), getColumn());
2327 if (!mIsSetReversible)
2328 {
2329 logError(AllowedAttributesOnReaction, level, version,
2330 "The required attribute 'reversible' is missing from the "
2331 + elplusid + ".");
2332 }
2333
2334 //
2335 // fast: boolean { use="required" } (L3v1 only)
2336 //
2337 if (version == 1)
2338 {
2339 mIsSetFast = attributes.readInto("fast", mFast, getErrorLog(),
2340 false, getLine(), getColumn());
2341 if (!mIsSetFast)
2342 {
2343 logError(AllowedAttributesOnReaction, level, version,
2344 "The required attribute 'fast' is missing from the "
2345 + elplusid + ".");
2346 }
2347 }
2348
2349 //
2350 // name: string { use="optional" } (L2v1 ->)
2351 //
2352 // for l3v2 sbase will read this
2353 if (version == 1)
2354 {
2355 attributes.readInto("name", mName, getErrorLog(), false,
2356 getLine(), getColumn());
2357 }
2358
2359 //
2360 // compartment: string { use="optional" } (L3v1 -> )
2361 //
2362 assigned = attributes.readInto("compartment", mCompartment, getErrorLog(), false, getLine(), getColumn());
2363 if (assigned && mCompartment.size() == 0)
2364 {
2365 logEmptyString("compartment", level, version, "<reaction>");
2366 }
2367 if (!SyntaxChecker::isValidInternalSId(mCompartment))
2368 logError(InvalidIdSyntax, getLevel(), getVersion(),
2369 "The " + elplusid + " has a 'compartment' with a value of '" + mCompartment
2370 + "' which does not conform .");
2371 }
2372 /** @endcond */
2373
2374
2375 /** @cond doxygenLibsbmlInternal */
2376 /*
2377 * Subclasses should override this method to write their XML attributes
2378 * to the XMLOutputStream. Be sure to call your parent's implementation
2379 * of this method as well.
2380 */
2381 void
writeAttributes(XMLOutputStream & stream) const2382 Reaction::writeAttributes (XMLOutputStream& stream) const
2383 {
2384 SBase::writeAttributes(stream);
2385
2386 const unsigned int level = getLevel ();
2387 const unsigned int version = getVersion();
2388
2389 //
2390 // sboTerm: SBOTerm { use="optional" } (L2v2->)
2391 //
2392 // sboTerm for L2V3 or later is written in SBase::writeAttributes()
2393 //
2394 if ( (level == 2) && (version == 2) )
2395 {
2396 SBO::writeTerm(stream, mSBOTerm);
2397 }
2398
2399 // for L3V2 and above SBase will write this out
2400 if (level < 3 || (level == 3 && version == 1))
2401 {
2402 //
2403 // name: SName { use="required" } (L1v1, L1v2)
2404 // id: SId { use="required" } (L2v1, L2v2)
2405 //
2406 const string id = (level == 1) ? "name" : "id";
2407 stream.writeAttribute(id, mId);
2408
2409 //
2410 // name: string { use="optional" } (L2v1->)
2411 //
2412 if (level > 1) stream.writeAttribute("name", mName);
2413 }
2414
2415 //
2416 // reversible: boolean { use="optional" default="true" }
2417 // (L1v1, L1v2, L2v1-> )
2418 // reversible: boolean { use="required"} (L3v1->)
2419 //
2420 if (level < 3)
2421 {
2422 if (mReversible != true || isExplicitlySetReversible())
2423 stream.writeAttribute("reversible", mReversible);
2424 }
2425 else
2426 {
2427 // in L3 only write it out if it has been set
2428 if (isSetReversible())
2429 stream.writeAttribute("reversible", mReversible);
2430 }
2431
2432 //
2433 // fast: boolean { use="optional" default="false" } (L1v1, L1v2, L2v1 ->)
2434 // fast: boolean { use="required" } (L3v1-> )
2435 //
2436 if (level < 3)
2437 {
2438 if (mIsSetFast)
2439 {
2440 if (isExplicitlySetFast() || level != 1 || mFast != false)
2441 stream.writeAttribute("fast", mFast);
2442 }
2443 }
2444 else
2445 {
2446 // in L3 only write it out if it has been set
2447 if (version == 1 && isSetFast())
2448 stream.writeAttribute("fast", mFast);
2449 }
2450
2451 if (level > 2)
2452 {
2453 stream.writeAttribute("compartment", mCompartment);
2454 }
2455
2456 //
2457 // (EXTENSION)
2458 //
2459 SBase::writeExtensionAttributes(stream);
2460 }
2461 /** @endcond */
2462
2463
2464 /** @cond doxygenLibsbmlInternal */
2465 /*
2466 * Subclasses should override this method to write out their contained
2467 * SBML objects as XML elements. Be sure to call your parent's
2468 * implementation of this method as well.
2469 */
2470 void
writeElements(XMLOutputStream & stream) const2471 Reaction::writeElements (XMLOutputStream& stream) const
2472 {
2473 SBase::writeElements(stream);
2474
2475 const unsigned int level = getLevel();
2476
2477 if (getLevel() == 3 && getVersion() > 1)
2478 {
2479 if (mReactants.hasOptionalElements() == true ||
2480 mReactants.hasOptionalAttributes() == true ||
2481 mReactants.isExplicitlyListed())
2482 {
2483 mReactants.write(stream);
2484 }
2485
2486 if (mProducts.hasOptionalElements() == true ||
2487 mProducts.hasOptionalAttributes() == true ||
2488 mProducts.isExplicitlyListed())
2489 {
2490 mProducts.write(stream);
2491 }
2492
2493 if (mModifiers.hasOptionalElements() == true ||
2494 mModifiers.hasOptionalAttributes() == true ||
2495 mModifiers.isExplicitlyListed())
2496 {
2497 mModifiers.write(stream);
2498 }
2499 }
2500 else
2501 {
2502 // use original code
2503 if (getNumReactants () > 0) mReactants.write(stream);
2504 if (getNumProducts () > 0) mProducts .write(stream);
2505
2506 if (level > 1 && getNumModifiers () > 0) mModifiers.write(stream);
2507 }
2508
2509 if (mKineticLaw != NULL) mKineticLaw->write(stream);
2510
2511 //
2512 // (EXTENSION)
2513 //
2514 SBase::writeExtensionElements(stream);
2515 }
2516 /** @endcond */
2517
2518
2519 /*
2520 * Creates a new ListOfReactions items.
2521 */
ListOfReactions(unsigned int level,unsigned int version)2522 ListOfReactions::ListOfReactions (unsigned int level, unsigned int version)
2523 : ListOf(level,version)
2524 {
2525 }
2526
2527
2528 /*
2529 * Creates a new ListOfReactions items.
2530 */
ListOfReactions(SBMLNamespaces * sbmlns)2531 ListOfReactions::ListOfReactions (SBMLNamespaces* sbmlns)
2532 : ListOf(sbmlns)
2533 {
2534 loadPlugins(sbmlns);
2535 }
2536
2537
2538 /*
2539 * @return a (deep) copy of this ListOfReactions.
2540 */
2541 ListOfReactions*
clone() const2542 ListOfReactions::clone () const
2543 {
2544 return new ListOfReactions(*this);
2545 }
2546
2547
2548 /*
2549 * @return the typecode (int) of SBML objects contained in this ListOf or
2550 * SBML_UNKNOWN (default).
2551 */
2552 int
getItemTypeCode() const2553 ListOfReactions::getItemTypeCode () const
2554 {
2555 return SBML_REACTION;
2556 }
2557
2558
2559 /*
2560 * @return the name of this element ie "listOfReactions".
2561 */
2562 const string&
getElementName() const2563 ListOfReactions::getElementName () const
2564 {
2565 static const string name = "listOfReactions";
2566 return name;
2567 }
2568
2569
2570 /* return nth item in list */
2571 Reaction *
get(unsigned int n)2572 ListOfReactions::get(unsigned int n)
2573 {
2574 return static_cast<Reaction*>(ListOf::get(n));
2575 }
2576
2577
2578 /* return nth item in list */
2579 const Reaction *
get(unsigned int n) const2580 ListOfReactions::get(unsigned int n) const
2581 {
2582 return static_cast<const Reaction*>(ListOf::get(n));
2583 }
2584
2585
2586 /* return item by id */
2587 Reaction*
get(const std::string & sid)2588 ListOfReactions::get (const std::string& sid)
2589 {
2590 return const_cast<Reaction*>(
2591 static_cast<const ListOfReactions&>(*this).get(sid) );
2592 }
2593
2594
2595 /* return item by id */
2596 const Reaction*
get(const std::string & sid) const2597 ListOfReactions::get (const std::string& sid) const
2598 {
2599 vector<SBase*>::const_iterator result;
2600
2601 result = find_if( mItems.begin(), mItems.end(), IdEq<Reaction>(sid) );
2602 return (result == mItems.end()) ? NULL : static_cast <Reaction*> (*result);
2603 }
2604
2605
2606 /* Removes the nth item from this list */
2607 Reaction*
remove(unsigned int n)2608 ListOfReactions::remove (unsigned int n)
2609 {
2610 return static_cast<Reaction*>(ListOf::remove(n));
2611 }
2612
2613
2614 /* Removes item in this list by id */
2615 Reaction*
remove(const std::string & sid)2616 ListOfReactions::remove (const std::string& sid)
2617 {
2618 SBase* item = NULL;
2619 vector<SBase*>::iterator result;
2620
2621 result = find_if( mItems.begin(), mItems.end(), IdEq<Reaction>(sid) );
2622
2623 if (result != mItems.end())
2624 {
2625 item = *result;
2626 mItems.erase(result);
2627 }
2628
2629
2630 return static_cast <Reaction*> (item);
2631 }
2632
2633
2634 /** @cond doxygenLibsbmlInternal */
2635 /*
2636 * @return the ordinal position of the element with respect to its siblings
2637 * or -1 (default) to indicate the position is not significant.
2638 */
2639 int
getElementPosition() const2640 ListOfReactions::getElementPosition () const
2641 {
2642 return 11;
2643 }
2644 /** @endcond */
2645
2646
2647 /** @cond doxygenLibsbmlInternal */
2648 /*
2649 * @return the SBML object corresponding to next XMLToken in the
2650 * XMLInputStream or @c NULL if the token was not recognized.
2651 */
2652 SBase*
createObject(XMLInputStream & stream)2653 ListOfReactions::createObject (XMLInputStream& stream)
2654 {
2655 const string& name = stream.peek().getName();
2656 SBase* object = NULL;
2657
2658 if (name == "reaction")
2659 {
2660 try
2661 {
2662 object = new Reaction(getSBMLNamespaces());
2663 }
2664 catch (SBMLConstructorException*)
2665 {
2666 object = new Reaction(SBMLDocument::getDefaultLevel(),
2667 SBMLDocument::getDefaultVersion());
2668 }
2669 catch ( ... )
2670 {
2671 object = new Reaction(SBMLDocument::getDefaultLevel(),
2672 SBMLDocument::getDefaultVersion());
2673 }
2674
2675 if (object != NULL) mItems.push_back(object);
2676 }
2677
2678 return object;
2679 }
2680 /** @endcond */
2681
2682 #endif /* __cplusplus */
2683
2684
2685 /** @cond doxygenIgnored */
2686 LIBSBML_EXTERN
2687 Reaction_t *
Reaction_create(unsigned int level,unsigned int version)2688 Reaction_create (unsigned int level, unsigned int version)
2689 {
2690 try
2691 {
2692 Reaction* obj = new Reaction(level,version);
2693 return obj;
2694 }
2695 catch (SBMLConstructorException)
2696 {
2697 return NULL;
2698 }
2699 }
2700
2701
2702 LIBSBML_EXTERN
2703 Reaction_t *
Reaction_createWithNS(SBMLNamespaces_t * sbmlns)2704 Reaction_createWithNS (SBMLNamespaces_t* sbmlns)
2705 {
2706 try
2707 {
2708 Reaction* obj = new Reaction(sbmlns);
2709 return obj;
2710 }
2711 catch (SBMLConstructorException)
2712 {
2713 return NULL;
2714 }
2715 }
2716
2717 LIBSBML_EXTERN
2718 void
Reaction_free(Reaction_t * r)2719 Reaction_free (Reaction_t *r)
2720 {
2721 delete r;
2722 }
2723
2724
2725 LIBSBML_EXTERN
2726 Reaction_t *
Reaction_clone(const Reaction_t * r)2727 Reaction_clone (const Reaction_t *r)
2728 {
2729 return (r != NULL) ? static_cast<Reaction*>( r->clone() ) : NULL;
2730 }
2731
2732
2733 LIBSBML_EXTERN
2734 void
Reaction_initDefaults(Reaction_t * r)2735 Reaction_initDefaults (Reaction_t *r)
2736 {
2737 if (r != NULL)
2738 r->initDefaults();
2739 }
2740
2741
2742 LIBSBML_EXTERN
2743 const XMLNamespaces_t *
Reaction_getNamespaces(Reaction_t * r)2744 Reaction_getNamespaces(Reaction_t *r)
2745 {
2746 return (r != NULL) ? r->getNamespaces() : NULL;
2747 }
2748
2749 LIBSBML_EXTERN
2750 const char *
Reaction_getId(const Reaction_t * r)2751 Reaction_getId (const Reaction_t *r)
2752 {
2753 return (r != NULL && r->isSetId()) ? r->getId().c_str() : NULL;
2754 }
2755
2756
2757 LIBSBML_EXTERN
2758 const char *
Reaction_getName(const Reaction_t * r)2759 Reaction_getName (const Reaction_t *r)
2760 {
2761 return (r != NULL && r->isSetName()) ? r->getName().c_str() : NULL;
2762 }
2763
2764
2765 LIBSBML_EXTERN
2766 KineticLaw_t *
Reaction_getKineticLaw(Reaction_t * r)2767 Reaction_getKineticLaw (Reaction_t *r)
2768 {
2769 return (r != NULL) ? r->getKineticLaw() : NULL;
2770 }
2771
2772
2773 LIBSBML_EXTERN
2774 int
Reaction_getReversible(const Reaction_t * r)2775 Reaction_getReversible (const Reaction_t *r)
2776 {
2777 return (r != NULL) ? static_cast<int>( r->getReversible() ) : 0;
2778 }
2779
2780
2781 LIBSBML_EXTERN
2782 int
Reaction_getFast(const Reaction_t * r)2783 Reaction_getFast (const Reaction_t *r)
2784 {
2785 return (r != NULL) ? static_cast<int>( r->getFast() ) : 0;
2786 }
2787
2788
2789 LIBSBML_EXTERN
2790 const char *
Reaction_getCompartment(const Reaction_t * r)2791 Reaction_getCompartment (const Reaction_t *r)
2792 {
2793 return (r != NULL && r->isSetCompartment()) ?
2794 r->getCompartment().c_str() : NULL;
2795 }
2796
2797
2798 LIBSBML_EXTERN
2799 int
Reaction_isSetId(const Reaction_t * r)2800 Reaction_isSetId (const Reaction_t *r)
2801 {
2802 return (r != NULL) ? static_cast<int>( r->isSetId() ) : 0;
2803 }
2804
2805
2806 LIBSBML_EXTERN
2807 int
Reaction_isSetName(const Reaction_t * r)2808 Reaction_isSetName (const Reaction_t *r)
2809 {
2810 return (r != NULL) ? static_cast<int>( r->isSetName() ) : 0;
2811 }
2812
2813
2814 LIBSBML_EXTERN
2815 int
Reaction_isSetKineticLaw(const Reaction_t * r)2816 Reaction_isSetKineticLaw (const Reaction_t *r)
2817 {
2818 return (r != NULL) ? static_cast<int>( r->isSetKineticLaw() ) : 0;
2819 }
2820
2821
2822 LIBSBML_EXTERN
2823 int
Reaction_isSetFast(const Reaction_t * r)2824 Reaction_isSetFast (const Reaction_t *r)
2825 {
2826 return (r != NULL) ? static_cast<int>( r->isSetFast() ) : 0;
2827 }
2828
2829
2830 LIBSBML_EXTERN
2831 int
Reaction_isSetCompartment(const Reaction_t * r)2832 Reaction_isSetCompartment (const Reaction_t *r)
2833 {
2834 return (r != NULL) ? static_cast<int>( r->isSetCompartment() ) : 0;
2835 }
2836
2837
2838 LIBSBML_EXTERN
2839 int
Reaction_isSetReversible(const Reaction_t * r)2840 Reaction_isSetReversible (const Reaction_t *r)
2841 {
2842 return (r != NULL) ? static_cast<int>( r->isSetReversible() ) : 0;
2843 }
2844
2845
2846 LIBSBML_EXTERN
2847 int
Reaction_setId(Reaction_t * r,const char * sid)2848 Reaction_setId (Reaction_t *r, const char *sid)
2849 {
2850 if (r != NULL)
2851 return (sid == NULL) ? r->setId("") : r->setId(sid);
2852 else
2853 return LIBSBML_INVALID_OBJECT;
2854 }
2855
2856
2857 LIBSBML_EXTERN
2858 int
Reaction_setName(Reaction_t * r,const char * name)2859 Reaction_setName (Reaction_t *r, const char *name)
2860 {
2861 if (r != NULL)
2862 return (name == NULL) ? r->unsetName() : r->setName(name);
2863 else
2864 return LIBSBML_INVALID_OBJECT;
2865 }
2866
2867
2868 LIBSBML_EXTERN
2869 int
Reaction_setKineticLaw(Reaction_t * r,const KineticLaw_t * kl)2870 Reaction_setKineticLaw (Reaction_t *r, const KineticLaw_t *kl)
2871 {
2872 if (r != NULL)
2873 return (kl == NULL) ? r->unsetKineticLaw() : r->setKineticLaw(kl);
2874 else
2875 return LIBSBML_INVALID_OBJECT;
2876 }
2877
2878
2879 LIBSBML_EXTERN
2880 int
Reaction_setReversible(Reaction_t * r,int value)2881 Reaction_setReversible (Reaction_t *r, int value)
2882 {
2883 if (r != NULL)
2884 return r->setReversible( static_cast<bool>(value) );
2885 else
2886 return LIBSBML_INVALID_OBJECT;
2887 }
2888
2889
2890 LIBSBML_EXTERN
2891 int
Reaction_setFast(Reaction_t * r,int value)2892 Reaction_setFast (Reaction_t *r, int value)
2893 {
2894 if (r != NULL)
2895 return r->setFast( static_cast<bool>(value) );
2896 else
2897 return LIBSBML_INVALID_OBJECT;
2898 }
2899
2900
2901 LIBSBML_EXTERN
2902 int
Reaction_setCompartment(Reaction_t * r,const char * compartment)2903 Reaction_setCompartment (Reaction_t *r, const char *compartment)
2904 {
2905 if (r != NULL)
2906 return (compartment == NULL) ? r->unsetCompartment() :
2907 r->setCompartment(compartment);
2908 else
2909 return LIBSBML_INVALID_OBJECT;
2910 }
2911
2912
2913 LIBSBML_EXTERN
2914 int
Reaction_unsetName(Reaction_t * r)2915 Reaction_unsetName (Reaction_t *r)
2916 {
2917 if (r != NULL)
2918 return r->unsetName();
2919 else
2920 return LIBSBML_INVALID_OBJECT;
2921 }
2922
2923
2924 LIBSBML_EXTERN
2925 int
Reaction_unsetCompartment(Reaction_t * r)2926 Reaction_unsetCompartment (Reaction_t *r)
2927 {
2928 if (r != NULL)
2929 return r->unsetCompartment();
2930 else
2931 return LIBSBML_INVALID_OBJECT;
2932 }
2933
2934
2935 LIBSBML_EXTERN
2936 int
Reaction_unsetKineticLaw(Reaction_t * r)2937 Reaction_unsetKineticLaw (Reaction_t *r)
2938 {
2939 if (r != NULL)
2940 return r->unsetKineticLaw();
2941 else
2942 return LIBSBML_INVALID_OBJECT;
2943 }
2944
2945
2946 LIBSBML_EXTERN
2947 int
Reaction_unsetFast(Reaction_t * r)2948 Reaction_unsetFast (Reaction_t *r)
2949 {
2950 if (r != NULL)
2951 return r->unsetFast();
2952 else
2953 return LIBSBML_INVALID_OBJECT;
2954 }
2955
2956
2957 LIBSBML_EXTERN
2958 int
Reaction_unsetReversible(Reaction_t * r)2959 Reaction_unsetReversible (Reaction_t *r)
2960 {
2961 if (r != NULL)
2962 return r->unsetReversible();
2963 else
2964 return LIBSBML_INVALID_OBJECT;
2965 }
2966
2967
2968 LIBSBML_EXTERN
2969 int
Reaction_hasRequiredAttributes(Reaction_t * r)2970 Reaction_hasRequiredAttributes(Reaction_t *r)
2971 {
2972 return (r != NULL) ? static_cast<int>(r->hasRequiredAttributes()) : 0;
2973 }
2974
2975
2976 LIBSBML_EXTERN
2977 int
Reaction_addReactant(Reaction_t * r,const SpeciesReference_t * sr)2978 Reaction_addReactant (Reaction_t *r, const SpeciesReference_t *sr)
2979 {
2980 if (r != NULL)
2981 return r->addReactant( static_cast<const SpeciesReference*>(sr) );
2982 else
2983 return LIBSBML_INVALID_OBJECT;
2984 }
2985
2986
2987 LIBSBML_EXTERN
2988 int
Reaction_addReactantBySpecies(Reaction_t * r,const Species_t * s,double stoichiometry,const char * id,int constant)2989 Reaction_addReactantBySpecies (Reaction_t *r, const Species_t *s,
2990 double stoichiometry, const char *id,
2991 int constant)
2992 {
2993 if (r != NULL)
2994 return r->addReactant( static_cast<const Species*>(s), stoichiometry, id,
2995 constant);
2996 else
2997 return LIBSBML_INVALID_OBJECT;
2998 }
2999
3000
3001 LIBSBML_EXTERN
3002 int
Reaction_addProduct(Reaction_t * r,const SpeciesReference_t * sr)3003 Reaction_addProduct (Reaction_t *r, const SpeciesReference_t *sr)
3004 {
3005 if (r != NULL)
3006 return r->addProduct( static_cast<const SpeciesReference*>(sr) );
3007 else
3008 return LIBSBML_INVALID_OBJECT;
3009 }
3010
3011
3012 LIBSBML_EXTERN
3013 int
Reaction_addProductBySpecies(Reaction_t * r,const Species_t * s,double stoichiometry,const char * id,int constant)3014 Reaction_addProductBySpecies (Reaction_t *r, const Species_t *s,
3015 double stoichiometry, const char *id,
3016 int constant)
3017 {
3018 if (r != NULL)
3019 return r->addProduct( static_cast<const Species*>(s), stoichiometry, id,
3020 constant);
3021 else
3022 return LIBSBML_INVALID_OBJECT;
3023 }
3024
3025
3026 LIBSBML_EXTERN
3027 int
Reaction_addModifier(Reaction_t * r,const SpeciesReference_t * msr)3028 Reaction_addModifier (Reaction_t *r, const SpeciesReference_t *msr)
3029 {
3030 if (r != NULL)
3031 {
3032 if (msr == NULL || msr->isModifier())
3033 {
3034 return r->addModifier(static_cast<const ModifierSpeciesReference*>(msr) );
3035 }
3036 else
3037 {
3038 return LIBSBML_INVALID_ATTRIBUTE_VALUE;
3039 }
3040 }
3041 else
3042 return LIBSBML_INVALID_OBJECT;
3043 }
3044
3045
3046 LIBSBML_EXTERN
3047 int
Reaction_addModifierBySpecies(Reaction_t * r,const Species_t * s,const char * id)3048 Reaction_addModifierBySpecies (Reaction_t *r, const Species_t *s,
3049 const char *id)
3050 {
3051 if (r != NULL)
3052 return r->addModifier( static_cast<const Species*>(s), id);
3053 else
3054 return LIBSBML_INVALID_OBJECT;
3055 }
3056
3057
3058 LIBSBML_EXTERN
3059 SpeciesReference_t *
Reaction_createReactant(Reaction_t * r)3060 Reaction_createReactant (Reaction_t *r)
3061 {
3062 return (r != NULL) ? r->createReactant() : NULL;
3063 }
3064
3065
3066 LIBSBML_EXTERN
3067 SpeciesReference_t *
Reaction_createProduct(Reaction_t * r)3068 Reaction_createProduct (Reaction_t *r)
3069 {
3070 return (r != NULL) ? r->createProduct() : NULL;
3071 }
3072
3073
3074 LIBSBML_EXTERN
3075 SpeciesReference_t *
Reaction_createModifier(Reaction_t * r)3076 Reaction_createModifier (Reaction_t *r)
3077 {
3078 return (r != NULL) ? r->createModifier() : NULL;
3079 }
3080
3081
3082 LIBSBML_EXTERN
3083 KineticLaw_t *
Reaction_createKineticLaw(Reaction_t * r)3084 Reaction_createKineticLaw (Reaction_t *r)
3085 {
3086 return (r != NULL) ? r->createKineticLaw() : NULL;
3087 }
3088
3089
3090 LIBSBML_EXTERN
3091 ListOf_t *
Reaction_getListOfReactants(Reaction_t * r)3092 Reaction_getListOfReactants (Reaction_t *r)
3093 {
3094 return (r != NULL) ? r->getListOfReactants() : NULL;
3095 }
3096
3097
3098 LIBSBML_EXTERN
3099 ListOf_t *
Reaction_getListOfProducts(Reaction_t * r)3100 Reaction_getListOfProducts (Reaction_t *r)
3101 {
3102 return (r != NULL) ? r->getListOfProducts() : NULL;
3103 }
3104
3105
3106 LIBSBML_EXTERN
3107 ListOf_t *
Reaction_getListOfModifiers(Reaction_t * r)3108 Reaction_getListOfModifiers (Reaction_t *r)
3109 {
3110 return (r != NULL) ? r->getListOfModifiers() : NULL;
3111 }
3112
3113
3114 LIBSBML_EXTERN
3115 SpeciesReference_t *
Reaction_getReactant(Reaction_t * r,unsigned int n)3116 Reaction_getReactant (Reaction_t *r, unsigned int n)
3117 {
3118 return (r != NULL) ? r->getReactant(n) : NULL;
3119 }
3120
3121
3122 LIBSBML_EXTERN
3123 SpeciesReference_t *
Reaction_getReactantBySpecies(Reaction_t * r,const char * species)3124 Reaction_getReactantBySpecies (Reaction_t *r, const char *species)
3125 {
3126 return (r != NULL && species != NULL) ? r->getReactant(species) : NULL;
3127 }
3128
3129
3130 LIBSBML_EXTERN
3131 SpeciesReference_t *
Reaction_getProduct(Reaction_t * r,unsigned int n)3132 Reaction_getProduct (Reaction_t *r, unsigned int n)
3133 {
3134 return (r != NULL) ? r->getProduct(n) : NULL;
3135 }
3136
3137
3138 LIBSBML_EXTERN
3139 SpeciesReference_t *
Reaction_getProductBySpecies(Reaction_t * r,const char * species)3140 Reaction_getProductBySpecies (Reaction_t *r, const char *species)
3141 {
3142 return (r != NULL && species != NULL) ? r->getProduct(species) : NULL;
3143 }
3144
3145
3146 LIBSBML_EXTERN
3147 SpeciesReference_t *
Reaction_getModifier(Reaction_t * r,unsigned int n)3148 Reaction_getModifier (Reaction_t *r, unsigned int n)
3149 {
3150 return (r != NULL) ? r->getModifier(n) : NULL;
3151 }
3152
3153
3154 LIBSBML_EXTERN
3155 SpeciesReference_t *
Reaction_getModifierBySpecies(Reaction_t * r,const char * species)3156 Reaction_getModifierBySpecies (Reaction_t *r, const char *species)
3157 {
3158 return (r != NULL && species != NULL) ? r->getModifier(species) : NULL;
3159 }
3160
3161 LIBSBML_EXTERN
3162 unsigned int
Reaction_getNumReactants(const Reaction_t * r)3163 Reaction_getNumReactants (const Reaction_t *r)
3164 {
3165 return (r != NULL) ? r->getNumReactants() : SBML_INT_MAX;
3166 }
3167
3168
3169 LIBSBML_EXTERN
3170 unsigned int
Reaction_getNumProducts(const Reaction_t * r)3171 Reaction_getNumProducts (const Reaction_t *r)
3172 {
3173 return (r != NULL) ? r->getNumProducts() : SBML_INT_MAX;
3174 }
3175
3176
3177 LIBSBML_EXTERN
3178 unsigned int
Reaction_getNumModifiers(const Reaction_t * r)3179 Reaction_getNumModifiers (const Reaction_t *r)
3180 {
3181 return (r != NULL) ? r->getNumModifiers() : SBML_INT_MAX;
3182 }
3183
3184
3185 LIBSBML_EXTERN
3186 SpeciesReference_t *
Reaction_removeReactant(Reaction_t * r,unsigned int n)3187 Reaction_removeReactant (Reaction_t *r, unsigned int n)
3188 {
3189 return (r != NULL) ? r->removeReactant(n) : NULL;
3190 }
3191
3192
3193 LIBSBML_EXTERN
3194 SpeciesReference_t *
Reaction_removeReactantBySpecies(Reaction_t * r,const char * species)3195 Reaction_removeReactantBySpecies (Reaction_t *r, const char *species)
3196 {
3197 if (r != NULL)
3198 return (species != NULL) ? r->removeReactant(species) : NULL;
3199 else
3200 return NULL;
3201 }
3202
3203
3204 LIBSBML_EXTERN
3205 SpeciesReference_t *
Reaction_removeProduct(Reaction_t * r,unsigned int n)3206 Reaction_removeProduct (Reaction_t *r, unsigned int n)
3207 {
3208 return (r != NULL) ? r->removeProduct(n) : NULL;
3209 }
3210
3211
3212 LIBSBML_EXTERN
3213 SpeciesReference_t *
Reaction_removeProductBySpecies(Reaction_t * r,const char * species)3214 Reaction_removeProductBySpecies (Reaction_t *r, const char *species)
3215 {
3216 if (r != NULL)
3217 return (species != NULL) ? r->removeProduct(species) : NULL;
3218 else
3219 return NULL;
3220 }
3221
3222
3223 LIBSBML_EXTERN
3224 SpeciesReference_t *
Reaction_removeModifier(Reaction_t * r,unsigned int n)3225 Reaction_removeModifier (Reaction_t *r, unsigned int n)
3226 {
3227 return (r != NULL) ? r->removeModifier(n) : NULL;
3228 }
3229
3230
3231 LIBSBML_EXTERN
3232 SpeciesReference_t *
Reaction_removeModifierBySpecies(Reaction_t * r,const char * species)3233 Reaction_removeModifierBySpecies (Reaction_t *r, const char *species)
3234 {
3235 if (r != NULL)
3236 return (species != NULL) ? r->removeModifier(species) : NULL;
3237 else
3238 return NULL;
3239 }
3240
3241
3242 LIBSBML_EXTERN
3243 Reaction_t *
ListOfReactions_getById(ListOf_t * lo,const char * sid)3244 ListOfReactions_getById (ListOf_t *lo, const char *sid)
3245 {
3246 if (lo != NULL)
3247 return (sid != NULL) ?
3248 static_cast <ListOfReactions *> (lo)->get(sid) : NULL;
3249 else
3250 return NULL;
3251 }
3252
3253
3254 LIBSBML_EXTERN
3255 Reaction_t *
ListOfReactions_removeById(ListOf_t * lo,const char * sid)3256 ListOfReactions_removeById (ListOf_t *lo, const char *sid)
3257 {
3258 if (lo != NULL)
3259 return (sid != NULL) ?
3260 static_cast <ListOfReactions *> (lo)->remove(sid) : NULL;
3261 else
3262 return NULL;
3263 }
3264 /** @endcond */
3265 LIBSBML_CPP_NAMESPACE_END
3266