File: regressionModel.C

package info (click to toggle)
ball 1.5.0%2Bgit20180813.37fc53c-6
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 239,888 kB
  • sloc: cpp: 326,149; ansic: 4,208; python: 2,303; yacc: 1,778; lex: 1,099; xml: 958; sh: 322; makefile: 95
file content (311 lines) | stat: -rw-r--r-- 8,428 bytes parent folder | download | duplicates (6)
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
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
// -*- Mode: C++; tab-width: 2; -*-
// vi: set ts=2:
//
// 

#include <BALL/QSAR/regressionModel.h>

using namespace std;

namespace BALL
{
	namespace QSAR
	{


		RegressionModel::RegressionModel(const QSARData& q) : Model(q) 
		{
			validation = new RegressionValidation(this);
			model_val = validation;
			offsets_.resize(0);
		}

		RegressionModel::~RegressionModel() 
		{
			delete validation;
		}

		void RegressionModel::operator = (const RegressionModel& m)
		{
			Model::operator = (m);
		}

		const Eigen::MatrixXd* RegressionModel::getTrainingResult() const
		{
			return &training_result_;
		}



		void RegressionModel::show()
		{
		// 	if (training_result_.rows() == 0)
		// 	{
		// 		throw Exception::InconsistentUsage(__FILE__, __LINE__, "Model must have been trained before the results can be saved displayed!"); 
		// 	}
			
			const Eigen::MatrixXd* coeffErrors = validation->getCoefficientStdErrors();
			bool sterr = 0;
			if (coeffErrors->cols() != 0)
			{
				sterr = 1;
			}
			cout<<type_<<"\t"<<data->getNoDescriptors()<<"\n";
			
			
			if (type_ == "ALL")
			{
				std::multiset<unsigned int>::iterator d_it = descriptor_IDs_.begin();
				for (int i = 0; d_it != descriptor_IDs_.end(); ++d_it, ++i)
				{
					cout<<*d_it<<"  "<<descriptor_names_[i]<<endl;
				}
			}
			else if (!descriptor_IDs_.empty())
			{
				std::multiset<unsigned int>::iterator d_it = descriptor_IDs_.begin();
				for (int i = 0; i < training_result_.rows(); i++)
				{
					if (type_ != "KPLS" && type_ != "KPCR" && type_ != "KPLS" && type_ != "GP")
					{
						cout<<String(*d_it)<<"\t"<<descriptor_names_[i]<<"\t";
						++d_it;
					}
					else
					{
						cout<<String(i)<<"    "<<substance_names_[i]<<"    ";
					}
					for (int j = 1; j <= training_result_.cols(); j++) 
					{
						cout<<training_result_(i+1, j)<<"\t";
					}
					for (int j = 1; j <= coeffErrors->cols(); j++)
					{
						cout<<(*coeffErrors)(i+1, j)<<"\t";
					}
					cout <<"\n";
				}
			}
			else
			{
				for (int i = 0; i < training_result_.rows(); i++)
				{
					if (type_ != "KPLS" && type_ != "KPCR" && type_ != "KPLS" && type_ != "GP")
					{
						cout<<String(i)<<"\t"<<descriptor_names_[i]<<"\t";
					}
					else
					{
						cout<<String(i)<<"\t"<<substance_names_[i]<<"\t";cout.flush();
					}
					for (int j = 1; j <= training_result_.cols(); j++) 
					{
						cout<<training_result_(i+1, j)<<"\t";
					}
					for (int j = 1; j <= coeffErrors->cols(); j++)
					{
						cout<<(*coeffErrors)(i+1, j)<<"\t";
					}
					cout <<"\n";
				}
			}

		}

		void RegressionModel::saveToFile(string filename)
		{
			if (data == 0)
			{
				cout<<"Error: no QSARData object assigned to model! Can not save model!"<<endl;
				return;
			}	
			
			bool trained = 1;
			if (training_result_.rows() == 0) trained = 0; 
			
			ofstream out(filename.c_str());
			
			const Eigen::MatrixXd* coeffErrors = validation->getCoefficientStdErrors();
			bool sterr = 0;
			if (coeffErrors->cols() != 0)
			{
				sterr = 1;
			}
			bool centered_data = 0;
			bool centered_y = 0;
			if (descriptor_transformations_.cols() != 0)
			{
				centered_data = 1;
				if (y_transformations_.cols() != 0)
				{
					centered_y = 1;
				}
			}
			
			int sel_features = descriptor_IDs_.size();
			if (sel_features == 0)
			{
				sel_features = data->getNoDescriptors();
			}
			
			int no_y = training_result_.cols();
			if (no_y == 0) no_y = y_transformations_.cols(); // correct no because transformation information will have to by read anyway when reading this model later ...
			
			out<<"# model-type_\tno of featues in input data\tselected featues\tno of response variables\tcentered descriptors?\tcentered response?\ttrained?"<<endl;
			out<<type_<<"\t"<<data->getNoDescriptors()<<"\t"<<sel_features<<"\t"<<no_y<<"\t"<<centered_data<<"\t"<<centered_y<<"\t"<<trained<<"\n\n";
			
			saveModelParametersToFile(out);
			saveResponseTransformationToFile(out); 
			RegressionModel::saveDescriptorInformationToFile(out); 
			out<<"# offsets"<<endl;
			out<<offsets_<<endl;
			
			out.close();
		}


		void RegressionModel::readFromFile(string filename)
		{
			ifstream input(filename.c_str()); 
			if (!input)
			{
				throw BALL::Exception::FileNotFound(__FILE__, __LINE__, filename);
			}	
			
			String line0;
			getline(input, line0);  // skip comment line 
			getline(input, line0);  // read read line containing model specification
			
			if (line0.getField(0, "\t") != type_)
			{	
				String e = "Wrong input data! Use training data file generated by a ";
				e = e + type_ + " model !";
				throw Exception::WrongDataType(__FILE__, __LINE__, e.c_str());
			}
			
			int no_descriptors = line0.getField(2, "\t").toInt();
			int no_y = line0.getField(3, "\t").toInt();
			bool centered_data = line0.getField(4, "\t").toInt();
			bool centered_y = line0.getField(5, "\t").toInt();
			bool trained = line0.getField(6, "\t").toInt();
			
			getline(input, line0);  // skip empty line
			readModelParametersFromFile(input);
			if (centered_y)
			{
				readResponseTransformationFromFile(input, no_y); 
			}
			readDescriptorInformationFromFile(input, no_descriptors, centered_data, no_y*trained); 
			getline(input, line0);  // skip empty line 
			getline(input, line0);  // skip comment line 
			if (input.eof()) offsets_.resize(0); 
			else readVector(offsets_, input, 1, no_y);
			
			input.close();
		}


		void RegressionModel::readDescriptorInformationFromFile(ifstream& input, int no_descriptors, bool transformation, int no_coefficients)
		{
			descriptor_names_.clear();
			if (transformation) descriptor_transformations_.resize(2, no_descriptors); 
			else descriptor_transformations_.resize(0, 0); 
			if (no_coefficients > 0) training_result_.resize(no_descriptors, no_coefficients); 
			else training_result_.resize(0, 0);
			String line;
			getline(input, line);  // skip comment line
			
			for (int i = 1; i <= no_descriptors; i++)
			{
				getline(input, line);
				unsigned int id = (unsigned int) line.getField(0, "\t").toInt();
				descriptor_IDs_.insert(id);
				descriptor_names_.push_back(line.getField(1, "\t"));
				int j = 2;
				for (; j < 2+no_coefficients; j++)
				{	
					training_result_(i, j-1) = line.getField(j, "\t").toDouble();
				}
				if (transformation)
				{
					descriptor_transformations_(1, i) = line.getField(j, "\t").toDouble(); 
					descriptor_transformations_(2, i) = line.getField(j+1, "\t").toDouble(); 
				}
			}
			getline(input, line);  // skip empty line	
		}

		void RegressionModel::saveDescriptorInformationToFile(ofstream& out)
		{
			out<<"# ID\tdescriptor-name\tcoefficient(s)\t";
			
			bool centered_data = (descriptor_transformations_.cols() > 0); 
			
			if (centered_data)
			{
				out<<"mean of desc.\tstddev of desc.\t";
			}
			if (stderr)
			{
				out<<"stderr(s) of coeff.";
			}
			out<<endl;
			
			const Eigen::MatrixXd* coeffErrors = validation->getCoefficientStdErrors();
			
			if (!descriptor_IDs_.empty())  // write descriptors and information about their transformation
			{
				descriptor_IDs_.begin();
				bool trained = (training_result_.rows() == descriptor_IDs_.size());
				
				std::multiset<unsigned int>::iterator d_it = descriptor_IDs_.begin();
				for (unsigned int i = 0; i < descriptor_IDs_.size(); i++, ++d_it)
				{
					out<<String(*d_it)<<"\t"<<descriptor_names_[i]<<"\t";
					if (trained)
					{
						for (int j = 1; j <= training_result_.cols(); j++) 
						{
							out<<training_result_(i+1, j)<<"\t";
						}
					}
					if (centered_data)
					{
						out<<descriptor_transformations_(1, i+1)<<"\t"<<descriptor_transformations_(2, i+1)<<"\t"; 
					}
					for (int j = 1; j <= coeffErrors->cols(); j++)
					{
						out<<(*coeffErrors)(i+1, j)<<"\t";
					}
					out <<"\n";
				}
			}
			else
			{
				bool trained = (training_result_.rows() == descriptor_names_.size());
				for (unsigned int i = 0; i < descriptor_names_.size(); i++)
				{
					out<<String(i)<<"\t"<<descriptor_names_[i]<<"\t";

					if (trained)
					{
						for (int j = 1; j <= training_result_.cols(); j++) 
						{
							out<<training_result_(i+1, j)<<"\t";
						}
					}
					if (centered_data)
					{
						out<<descriptor_transformations_(1, i+1)<<"\t"<<descriptor_transformations_(2, i+1)<<"\t"; 
					}
					for (int j = 1; j <= coeffErrors->cols(); j++)
					{
						out<<(*coeffErrors)(i+1, j)<<"\t";
					}
					out <<"\n";
				}
			}
			out<<endl;
		}
	}
}