File: isValidSBML_Model.m

package info (click to toggle)
sbmltoolbox 4.1.0-5.1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 7,032 kB
  • sloc: xml: 2,438; makefile: 8; sh: 7
file content (242 lines) | stat: -rw-r--r-- 6,930 bytes parent folder | download | duplicates (5)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
function [valid, message] = isValidSBML_Model(SBMLStructure)
% [valid, message] = isValidSBML_Model(SBMLModel)
%
% Takes
%
% 1. SBMLModel, an SBML Model structure
%
% Returns
%
% 1. valid = 
%   - 1, if the structure represents
%        a MATLAB_SBML Model structure of the appropriate
%        level and version
%   - 0, otherwise
% 2. a message explaining any failure
%
% *NOTE:* The fields present in a MATLAB_SBML Model structure of the appropriate
% level and version can be found using getModelFieldnames(level, version)

%<!---------------------------------------------------------------------------
% This file is part of SBMLToolbox.  Please visit http://sbml.org for more
% information about SBML, and the latest version of SBMLToolbox.
%
% Copyright (C) 2009-2012 jointly by the following organizations: 
%     1. California Institute of Technology, Pasadena, CA, USA
%     2. EMBL European Bioinformatics Institute (EBML-EBI), Hinxton, UK
%
% Copyright (C) 2006-2008 jointly by the following organizations: 
%     1. California Institute of Technology, Pasadena, CA, USA
%     2. University of Hertfordshire, Hatfield, UK
%
% Copyright (C) 2003-2005 jointly by the following organizations: 
%     1. California Institute of Technology, Pasadena, CA, USA 
%     2. Japan Science and Technology Agency, Japan
%     3. University of Hertfordshire, Hatfield, UK
%
% SBMLToolbox is free software; you can redistribute it and/or modify it
% under the terms of the GNU Lesser General Public License as published by
% the Free Software Foundation.  A copy of the license agreement is provided
% in the file named "LICENSE.txt" included with this software distribution.
%----------------------------------------------------------------------- -->




%check the input arguments are appropriate

if (length(SBMLStructure) > 1)
  valid = 0;
  message = 'cannot deal with arrays of structures';
  return;
end;

if ~isempty(SBMLStructure)
  if isfield(SBMLStructure, 'SBML_level')
    level = SBMLStructure.SBML_level;
  else
    level = 3;
  end;
  if isfield(SBMLStructure, 'SBML_version')
    version = SBMLStructure.SBML_version;
  else
    version = 1;
  end;
else
  level = 3;
  version = 1;
end;

isValidLevelVersionCombination(level, version);

message = '';

% check that argument is a structure
valid = isstruct(SBMLStructure);

% check the typecode
typecode = 'SBML_MODEL';
if (valid == 1 && ~isempty(SBMLStructure))
  if isfield(SBMLStructure, 'typecode')
    if (strcmp(typecode, SBMLStructure.typecode) ~= 1)
      valid = 0;
      message = 'typecode mismatch';
      return;
    end;
  else
    valid = 0;
    message = 'missing typecode field';
    return;
  end;
end;


% check that structure contains all the necessary fields
[SBMLfieldnames, numFields] = getFieldnames('SBML_MODEL', level, version);

if (numFields ==0)
	valid = 0;
	message = 'invalid level/version';
end;

index = 1;
while (valid == 1 && index <= numFields)
	valid = isfield(SBMLStructure, char(SBMLfieldnames(index)));
	if (valid == 0);
		message = sprintf('%s field missing', char(SBMLfieldnames(index)));
	end;
	index = index + 1;
end;

%check that any nested structures are appropriate

% functionDefinitions
if (valid == 1 && level > 1)
  index = 1;
  while (valid == 1 && index <= length(SBMLStructure.functionDefinition))
    [valid, message] = isSBML_FunctionDefinition( ...
                                  SBMLStructure.functionDefinition(index), ...
                                  level, version);
    index = index + 1;
  end;
end;

% unitDefinitions
if (valid == 1)
  index = 1;
  while (valid == 1 && index <= length(SBMLStructure.unitDefinition))
    [valid, message] = isSBML_UnitDefinition( ...
                                  SBMLStructure.unitDefinition(index), ...
                                  level, version);
    index = index + 1;
  end;
end;

% compartments
if (valid == 1)
  index = 1;
  while (valid == 1 && index <= length(SBMLStructure.compartment))
    [valid, message] = isSBML_Compartment(SBMLStructure.compartment(index), ...
                                  level, version);
    index = index + 1;
  end;
end;

% species
if (valid == 1)
  index = 1;
  while (valid == 1 && index <= length(SBMLStructure.species))
    [valid, message] = isSBML_Species(SBMLStructure.species(index), ...
                                  level, version);
    index = index + 1;
  end;
end;

% compartmentTypes
if (valid == 1 && level == 2 && version > 1)
  index = 1;
  while (valid == 1 && index <= length(SBMLStructure.compartmentType))
    [valid, message] = isSBML_CompartmentType(SBMLStructure.compartmentType(index), ...
                                  level, version);
    index = index + 1;
  end;
end;

% speciesTypes
if (valid == 1 && level == 2 && version > 1)
  index = 1;
  while (valid == 1 && index <= length(SBMLStructure.speciesType))
    [valid, message] = isSBML_SpeciesType(SBMLStructure.speciesType(index), ...
                                  level, version);
    index = index + 1;
  end;
end;

% parameter
if (valid == 1)
  index = 1;
  while (valid == 1 && index <= length(SBMLStructure.parameter))
    [valid, message] = isSBML_Parameter(SBMLStructure.parameter(index), ...
                                  level, version);
    index = index + 1;
  end;
end;

% initialAssignment
if (valid == 1 && (level > 2 || (level == 2 && version > 1)))
  index = 1;
  while (valid == 1 && index <= length(SBMLStructure.initialAssignment))
    [valid, message] = isSBML_InitialAssignment( ...
                                  SBMLStructure.initialAssignment(index), ...
                                  level, version);
    index = index + 1;
  end;
end;

% rule
if (valid == 1)
  index = 1;
  while (valid == 1 && index <= length(SBMLStructure.rule))
    [valid, message] = isSBML_Rule(SBMLStructure.rule(index), ...
                                  level, version);
    index = index + 1;
  end;
end;

% constraints
if (valid == 1 && (level > 2 || (level == 2 && version > 1)))
  index = 1;
  while (valid == 1 && index <= length(SBMLStructure.constraint))
    [valid, message] = isSBML_Constraint( ...
                                  SBMLStructure.constraint(index), ...
                                  level, version);
    index = index + 1;
  end;
end;

% reaction
if (valid == 1)
  index = 1;
  while (valid == 1 && index <= length(SBMLStructure.reaction))
    [valid, message] = isSBML_Reaction(SBMLStructure.reaction(index), ...
                                  level, version);
    index = index + 1;
  end;
end;

% event
if (valid == 1 && level > 1)
  index = 1;
  while (valid == 1 && index <= length(SBMLStructure.event))
    [valid, message] = isSBML_Event(SBMLStructure.event(index), ...
                                  level, version);
    index = index + 1;
  end;
end;


% report failure
if (valid == 0)
	message = sprintf('Invalid Model structure\n%s\n', message);
end;