Skip to content
Snippets Groups Projects
Commit 440f385c authored by ardalan's avatar ardalan
Browse files

new Parameter for export

parent bcb360a9
No related branches found
No related tags found
No related merge requests found
......@@ -331,8 +331,7 @@ deformatGene<-function(idstr) {
return(idstr)
}
exportSBML<- function(morg=NULL,level=2,version=4,FbcLevel=0,filename="export.xml",validation=TRUE){
#morg<-modelorg[["modelorg"]]
exportSBML<- function(morg=NULL,level=2,version=4,FbcLevel=0,filename="export.xml",recoverExtMet=TRUE,printNotes=TRUE,printAnnos=TRUE,validation=FALSE ){
if(class(morg)!="modelorg"){
stop("morg has to be of class modelorg\n")
}
......@@ -382,12 +381,14 @@ exportSBML<- function(morg=NULL,level=2,version=4,FbcLevel=0,filename="export.xm
##EXCHANGE REACTIONS##
ex <- findExchReact(morg)
ex_react<-NULL
if(!is.null(ex))
# if recoverExtMet== FALSE => null for ex_react
if( (!is.null(ex)) && (recoverExtMet) )
{
if(!(all(diag(S(morg)[met_pos(ex), react_pos(ex)])==-1)))
stop("exchange reactions with Scoeff different than -1\n")
ex_react<-as.integer(react_pos(ex))
}
### Build wrapper for C Function #####
......@@ -415,27 +416,24 @@ exportSBML<- function(morg=NULL,level=2,version=4,FbcLevel=0,filename="export.xm
newsubS<- NULL
# TODO ATRRIBUTE NUR LESEN BEI NEUER SYBIL VERSION
if( .hasSlot(morg,"mod_attr") && .hasSlot(morg,"comp_attr") && .hasSlot(morg,"met_attr") && .hasSlot(morg,"react_attr") )
newSybil<-TRUE
else newSybil<-FALSE
### Start newSybil attr
if(newSybil)
{
if("notes" %in% colnames(mod_attr(morg))) mod_notes<-as.character(mod_attr(morg)[['notes']])
if("annotation" %in% colnames(mod_attr(morg))) mod_annotation<-as.character(mod_attr(morg)[['annotation']])
if(("notes" %in% colnames(mod_attr(morg))) && (printNotes) ) mod_notes<-as.character(mod_attr(morg)[['notes']])
if(("annotation" %in% colnames(mod_attr(morg))) && (printAnnos) ) mod_annotation<-as.character(mod_attr(morg)[['annotation']])
if("notes" %in% colnames(comp_attr(morg))) com_notes<-as.character(as.list((comp_attr(morg)[['notes']])))
if("annotation" %in% colnames(comp_attr(morg))) com_annotation<-as.character(as.list((comp_attr(morg)[['annotation']])))
if(("notes" %in% colnames(comp_attr(morg))) && (printNotes) ) com_notes<-as.character(as.list((comp_attr(morg)[['notes']])))
if(("annotation" %in% colnames(comp_attr(morg))) && (printAnnos) ) com_annotation<-as.character(as.list((comp_attr(morg)[['annotation']])))
if("charge" %in% colnames(met_attr(morg))) met_charge<- as.integer(as.list((met_attr(morg)[['charge']])))
if("chemicalFormula" %in% colnames(met_attr(morg))) met_formula<-as.character(as.list((met_attr(morg)[['chemicalFormula']])))
if("annotation" %in% colnames(met_attr(morg))) met_anno<-as.character(as.list((met_attr(morg)[['annotation']])))
if(("annotation" %in% colnames(met_attr(morg))) && (printAnnos)) met_anno<-as.character(as.list((met_attr(morg)[['annotation']])))
if("boundaryCondition" %in% colnames(met_attr(morg))) met_bnd<-as.logical(as.list((met_attr(morg)[['boundaryCondition']])))
if("notes" %in% colnames(met_attr(morg)))
if(("notes" %in% colnames(met_attr(morg))) && (printNotes) )
{ # delete Formular and charge from notes to do
met_notes<-as.character(as.list((met_attr(morg)[['notes']])))
if (!is.null(met_charge) || !is.null(met_formula))
......@@ -481,10 +479,10 @@ exportSBML<- function(morg=NULL,level=2,version=4,FbcLevel=0,filename="export.xm
}
}
if("annotation" %in% colnames(react_attr(morg))) react_anno<-as.character(as.list((react_attr(morg)[['annotation']])))
if(("annotation" %in% colnames(react_attr(morg))) && (printAnnos)) react_anno<-as.character(as.list((react_attr(morg)[['annotation']])))
# Merge Notes with "our" Notes and make sure gpr Rules from gpr
if("notes" %in% colnames(react_attr(morg)))
if(("notes" %in% colnames(react_attr(morg))) && (printNotes))
{
react_notes<-as.character(as.list((react_attr(morg)[['notes']])))
# using
......@@ -565,6 +563,7 @@ exportSBML<- function(morg=NULL,level=2,version=4,FbcLevel=0,filename="export.xm
react_notes,
react_anno,
ex_react,
as.logical(validation),
as.character(deformatGene(allgenes))
)
return (success)
......
......@@ -54,7 +54,7 @@ static const R_CallMethodDef callMethods[] = {
{"getSBMLCompartList", (DL_FUNC) &getSBMLCompartList, 1},
{"getSBMLSpeciesList", (DL_FUNC) &getSBMLSpeciesList, 1},
{"getSBMLReactionsList", (DL_FUNC) &getSBMLReactionsList, 1},
{"exportSBML", (DL_FUNC) &exportSBML, 33},
{"exportSBML", (DL_FUNC) &exportSBML, 34},
{"getSBMLFbcversion", (DL_FUNC) &getSBMLFbcversion, 1},
{NULL, NULL, 0}
};
......
......@@ -1605,7 +1605,7 @@ void ParseModtoAnno (SBase_t* comp , char* Mannocopy)
SEXP exportSBML (SEXP version, SEXP level,SEXP FbcLevel, SEXP filename,SEXP sybil_max, SEXP mod_desc, SEXP mod_name, SEXP mod_compart, SEXP met_id, SEXP met_name, SEXP met_comp, SEXP met_form,SEXP met_charge, SEXP react_id, SEXP react_name, SEXP react_rev, SEXP lowbnd, SEXP uppbnd, SEXP obj_coef, SEXP subSys, SEXP gpr, SEXP SMatrix, SEXP mod_notes, SEXP mod_anno, SEXP com_notes , SEXP com_anno, SEXP met_notes, SEXP met_anno, SEXP met_bnd , SEXP react_notes, SEXP react_anno, SEXP ex_react, SEXP allgenes)
SEXP exportSBML (SEXP version, SEXP level,SEXP FbcLevel, SEXP filename,SEXP sybil_max, SEXP mod_desc, SEXP mod_name, SEXP mod_compart, SEXP met_id, SEXP met_name, SEXP met_comp, SEXP met_form,SEXP met_charge, SEXP react_id, SEXP react_name, SEXP react_rev, SEXP lowbnd, SEXP uppbnd, SEXP obj_coef, SEXP subSys, SEXP gpr, SEXP SMatrix, SEXP mod_notes, SEXP mod_anno, SEXP com_notes , SEXP com_anno, SEXP met_notes, SEXP met_anno, SEXP met_bnd , SEXP react_notes, SEXP react_anno, SEXP ex_react,SEXP validation, SEXP allgenes)
{
//Varaibles from R
const char* fname = CHAR(STRING_ELT(filename, 0));
......@@ -1791,7 +1791,6 @@ SEXP exportSBML (SEXP version, SEXP level,SEXP FbcLevel, SEXP filename,SEXP sybi
Compartment_setId(comp,sName);
Compartment_setConstant(comp,1);
if( strcmp(sName,"BOUNDARY")==0 || strcmp(sName,"Boundary")==0 || strcmp(sName,"boundary")==0)hasBoundary=1;
if (!Rf_isNull(com_notes) && Rf_length(com_notes) > 1)
{
char *Cnotes = (char*) CHAR(STRING_ELT(com_notes, i));
......@@ -1819,7 +1818,7 @@ SEXP exportSBML (SEXP version, SEXP level,SEXP FbcLevel, SEXP filename,SEXP sybi
/* Boundary Compartment */
if(hasBoundary==0 && Rf_isNull(met_bnd) && Rf_length(met_bnd) <= 1 )
if(hasBoundary==0 && Rf_isNull(met_bnd) && Rf_length(met_bnd) <= 1 && !Rf_isNull(ex_react))
{
comp = Model_createCompartment(model);
Compartment_setId(comp,"BOUNDARY");
......@@ -1903,13 +1902,8 @@ SEXP exportSBML (SEXP version, SEXP level,SEXP FbcLevel, SEXP filename,SEXP sybi
if((Manno != NULL) && (Manno[0] != '\0' ))
{
//XMLNode_t * xmlannotation= RDFAnnotationParser_createAnnotation();
//XMLNode_t *rdf=RDFAnnotationParser_createRDFAnnotation(2,4);
//XMLNode_addChild(xmlannotation,(const XMLNode_t*) rdf);
//char *Mmetaid = (char*) CHAR(STRING_ELT(met_metaid, i));
SBase_setMetaId((SBase_t*)sp, CHAR(STRING_ELT(met_id, i)));
//SBase_appendAnnotation((SBase_t*)sp, xmlannotation );
// COPY STRING
char *Manno = (char*) CHAR(STRING_ELT(met_anno, i));
......@@ -2094,7 +2088,7 @@ SEXP exportSBML (SEXP version, SEXP level,SEXP FbcLevel, SEXP filename,SEXP sybi
SpeciesReference_setStoichiometry(spr, fabs(REAL(SMatrix)[hash]));
//is Exchange Reaction
if(isexchange==1)
if(isexchange==1 && !Rf_isNull(ex_react))
{
/* Create boundary Species */
sp = Model_createSpecies(model);
......@@ -2281,6 +2275,20 @@ SEXP exportSBML (SEXP version, SEXP level,SEXP FbcLevel, SEXP filename,SEXP sybi
}
}
// if VALIDATION TRUE FIND ERRORS
if(LOGICAL(validation)[0]==1)
{
unsigned int errors = 0;
errors = SBMLDocument_getNumErrors(sbmlDoc);
// errors += SBMLDocument_checkConsistency(sbmlDoc); warnings
if (errors > 0){
SBMLDocument_printErrors(sbmlDoc, stdout);
printf("\n");
}
else printf("No errors \n");
}
// write SBML file
int result = writeSBML(sbmlDoc, fname);
SEXP out = R_NilValue;
if (result) out = Rf_mkString(append_strings("Wrote file",fname," "));
......
......@@ -109,5 +109,5 @@ SEXP getSBMLSpeciesList(SEXP sbmlmod);
SEXP getSBMLReactionsList(SEXP sbmlmod);
/* export Modelorg to SBML*/
SEXP exportSBML (SEXP version, SEXP level,SEXP FbcLevel, SEXP filename,SEXP sybil_max, SEXP mod_desc, SEXP mod_name, SEXP mod_compart, SEXP met_id, SEXP met_name, SEXP met_comp, SEXP met_form,SEXP met_charge, SEXP react_id, SEXP react_name, SEXP react_rev, SEXP lowbnd, SEXP uppbnd, SEXP obj_coef, SEXP subSys, SEXP gpr, SEXP SMatrix, SEXP mod_notes, SEXP mod_anno, SEXP com_notes , SEXP com_anno, SEXP met_notes, SEXP met_anno, SEXP met_bnd , SEXP react_notes, SEXP react_anno, SEXP ex_react, SEXP allgenes);
SEXP exportSBML (SEXP version, SEXP level,SEXP FbcLevel, SEXP filename,SEXP sybil_max, SEXP mod_desc, SEXP mod_name, SEXP mod_compart, SEXP met_id, SEXP met_name, SEXP met_comp, SEXP met_form,SEXP met_charge, SEXP react_id, SEXP react_name, SEXP react_rev, SEXP lowbnd, SEXP uppbnd, SEXP obj_coef, SEXP subSys, SEXP gpr, SEXP SMatrix, SEXP mod_notes, SEXP mod_anno, SEXP com_notes , SEXP com_anno, SEXP met_notes, SEXP met_anno, SEXP met_bnd , SEXP react_notes, SEXP react_anno, SEXP ex_react,SEXP validation, SEXP allgenes);
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment