source: trunk/external-alma/components/SpectralComponents/SpectralList.cc@ 2987

Last change on this file since 2987 was 2980, checked in by Malte Marquarding, 10 years ago

Add a copy of casacore/components/SpectralComponents to external-alma directory to prepare for its removal from casacore-trunk. DOn;t activate in SConscript yet.

File size: 7.5 KB
Line 
1//# SpectralList.cc: A set of SpectralElements
2//# Copyright (C) 2001,2002
3//# Associated Universities, Inc. Washington DC, USA.
4//#
5//# This library is free software; you can redistribute it and/or modify it
6//# under the terms of the GNU Library General Public License as published by
7//# the Free Software Foundation; either version 2 of the License, or (at your
8//# option) any later version.
9//#
10//# This library is distributed in the hope that it will be useful, but WITHOUT
11//# ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
12//# FITNESS FOR A PARTICULAR PURPOSE. See the GNU Library General Public
13//# License for more details.
14//#
15//# You should have received a copy of the GNU Library General Public License
16//# along with this library; if not, write to the Free Software Foundation,
17//# Inc., 675 Massachusetts Ave, Cambridge, MA 02139, USA.
18//#
19//# Correspondence concerning AIPS++ should be addressed as follows:
20//# Internet email: aips2-request@nrao.edu.
21//# Postal address: AIPS++ Project Office
22//# National Radio Astronomy Observatory
23//# 520 Edgemont Road
24//# Charlottesville, VA 22903-2475 USA
25//#
26//# $Id: SpectralList.cc 21465 2014-06-19 05:56:56Z gervandiepen $
27
28//# Includes
29#include <components/SpectralComponents/SpectralList.h>
30
31#include <casa/Exceptions/Error.h>
32#include <casa/Containers/RecordInterface.h>
33#include <casa/Containers/Record.h>
34#include <components/SpectralComponents/GaussianSpectralElement.h>
35#include <casa/Utilities/PtrHolder.h>
36#include <components/SpectralComponents/SpectralElementFactory.h>
37
38#include <casa/iostream.h>
39
40namespace casa { //# NAMESPACE CASA - BEGIN
41
42//# Constructors
43SpectralList::SpectralList() :
44 nmax_p(0), list_p(0) {}
45
46SpectralList::SpectralList(uInt nmax) :
47 nmax_p(nmax), list_p(0) {
48}
49
50SpectralList::SpectralList(const SpectralElement &in) :
51 nmax_p(0), list_p(1) {
52 list_p[0] = in.clone();
53}
54
55SpectralList::SpectralList(const SpectralList &other) :
56 nmax_p(other.nmax_p), list_p(other.list_p.nelements()) {
57 for (uInt i=0; i<list_p.nelements(); i++) {
58 list_p[i] = other.list_p[i]->clone();
59 }
60}
61
62SpectralList::~SpectralList() {
63 clear();
64}
65
66SpectralList &SpectralList::operator=(const SpectralList &other) {
67 if (this != &other) {
68 clear();
69 nmax_p = other.nmax_p;
70 list_p.resize(other.list_p.nelements());
71 for (uInt i=0; i<list_p.nelements(); i++) {
72 list_p[i] = other.list_p[i]->clone();
73 }
74 }
75 return *this;
76}
77
78Double SpectralList::operator()(const Double x) const {
79 Double s(0);
80 for (uInt i=0; i<list_p.nelements(); i++) s += (*list_p[i])(x);
81 return s;
82}
83
84const SpectralElement* SpectralList::operator[](const uInt n) const {
85 if (n >= list_p.nelements()) {
86 throw(AipsError("SpectralList: Illegal index for element"));
87 }
88 return list_p[n];
89}
90
91SpectralElement* SpectralList::operator[](const uInt n) {
92 if (n >= list_p.nelements()) {
93 throw(AipsError("SpectralList: Illegal index for element"));
94 }
95 return list_p[n];
96}
97
98Bool SpectralList::add(const SpectralElement &in) {
99 uInt i = list_p.nelements();
100 if (nmax_p != 0 && i >= nmax_p) return False;
101 list_p.resize(i+1);
102 list_p[i] = in.clone();
103 return True;
104}
105
106Bool SpectralList::add(const SpectralList &in) {
107 for (uInt i=0; i<in.nelements(); i++) {
108 if (! add(*in[i])) {
109 return False;
110 }
111 }
112 return True;
113}
114
115void SpectralList::insert(const SpectralElement &in) {
116 uInt n = list_p.nelements();
117 uInt i;
118 for (i=0; i<n; i++) {
119 if (compar(in, *list_p[i]) > 0) {
120 break;
121 }
122 }
123 if (i == n) add(in);
124 else {
125 if (nmax_p != 0 && n >= nmax_p) {
126 delete list_p[n-1]; list_p[n-1] = 0;
127 } else {
128 list_p.resize(n+1);
129 list_p[n++] = 0;
130 }
131 for (uInt j=n-1; j>i; j--) list_p[j] = list_p[j-1];
132 if (in.getType() == SpectralElement::GAUSSIAN) {
133 const GaussianSpectralElement *gIn = dynamic_cast<const GaussianSpectralElement *>(&in);
134 list_p[i] = new GaussianSpectralElement(*gIn);
135 }
136 else {
137 // FIXME for other subclasses
138 list_p[i] = in.clone();
139 }
140 }
141}
142
143void SpectralList::insert(const SpectralList &in) {
144 for (uInt i=0; i<in.nelements(); i++) {
145 insert(*in[i]);
146 }
147}
148
149Bool SpectralList::set(const SpectralElement &in, const uInt which) {
150 uInt i = list_p.nelements();
151 if (nmax_p != 0 && which >= nmax_p) return False;
152 if (which > i) return False;
153 if (which == i) add(in);
154 delete list_p[which]; list_p[which] = 0;
155 list_p[which] = in.clone();
156 return True;
157}
158
159void SpectralList::clear() {
160 for (uInt i=0; i<list_p.nelements(); i++) {
161 delete list_p[i]; list_p[i] = 0;
162 }
163 list_p.resize(0, True);
164}
165
166void SpectralList::set(const uInt nmax) {
167 if (nmax != 0 && nmax < list_p.nelements()) {
168 for (uInt i=nmax; i<list_p.nelements(); i++) {
169 delete list_p[i]; list_p[i] = 0;
170 }
171 list_p.resize(nmax, True);
172 }
173 nmax_p = nmax;
174}
175
176Bool SpectralList::fromRecord (String& errMsg, const RecordInterface& container)
177{
178 clear();
179 for (uInt i=0; i<container.nfields(); i++) {
180 if (container.dataType(i)==TpRecord) {
181 const RecordInterface& rec = container.asRecord(i);
182 PtrHolder<SpectralElement> specEl(SpectralElementFactory::fromRecord(rec));
183 add(*specEl);
184 } else {
185 errMsg = String("Illegal record structure");
186 return False;
187 }
188 }
189 return True;
190}
191
192Bool SpectralList::toRecord(RecordInterface& container) const {
193 String errMsg;
194 for (uInt i=0; i<list_p.nelements(); i++) {
195 Record elRec;
196 list_p[i]->toRecord(elRec);
197 container.defineRecord(i, elRec);
198 }
199 return True;
200}
201
202
203void SpectralList::sort() {
204 uInt n = list_p.nelements();
205 if (n < 2) return;
206 SpectralElement *x;
207 for (uInt i=0; i<n-1; i++) {
208 for (uInt j=n-1; j>i; j--) {
209 if (compar(*list_p[j-1], *list_p[j]) < 0) {
210 x = list_p[j-1];
211 list_p[j-1] = list_p[j];
212 list_p[j] = x;
213 }
214 }
215 }
216}
217
218Int SpectralList::compar(
219 const SpectralElement &p1,
220 const SpectralElement &p2
221) const {
222 SpectralElement::Types p1Type = p1.getType();
223 SpectralElement::Types p2Type = p2.getType();
224 Double p1Amp = 0;
225 Double p2Amp = 0;
226 if (p1Type == SpectralElement::GAUSSIAN) {
227 const GaussianSpectralElement *g1 = dynamic_cast<const GaussianSpectralElement *>(&p1);
228 p1Amp = g1->getAmpl();
229 }
230 if (p2Type == SpectralElement::GAUSSIAN) {
231 const GaussianSpectralElement *g2 = dynamic_cast<const GaussianSpectralElement *>(&p2);
232 p2Amp = g2->getAmpl();
233 }
234 if (p1Amp > p2Amp) {
235 return 1;
236 }
237 else if (p1Amp < p2Amp) {
238 return -1;
239 }
240 else {
241 return 0;
242 }
243}
244
245ostream &operator<<(ostream &os, const SpectralList &lst) {
246 os << lst.nelements() << " in SpectralList:" << endl;
247 for (uInt i=0; i<lst.nelements(); i++) os << *lst[i];
248
249 return os;
250}
251
252} //# NAMESPACE CASA - END
253
254
255//# Cater for Double and Float
256#ifdef AIPS_NO_TEMPLATE_SRC
257#include <components/SpectralComponents/SpectralList2.tcc>
258
259namespace casa { //# NAMESPACE CASA - BEGIN
260template void SpectralList::residual<Double>(Vector<Double> &) const;
261template void SpectralList::residual<Float>(Vector<Float> &) const;
262template void SpectralList::evaluate<Double>(Vector<Double> &) const;
263template void SpectralList::evaluate<Float>(Vector<Float> &) const;
264template void SpectralList::residual<Double>(Vector<Double> &,
265 Vector<Double> const &) const;
266template void SpectralList::residual<Float>(Vector<Float> &,
267 Vector<Float> const &) const;
268template void SpectralList::evaluate<Double>(Vector<Double> &,
269 Vector<Double> const &) const;
270template void SpectralList::evaluate<Float>(Vector<Float> &,
271 Vector<Float> const &) const;
272} //# NAMESPACE CASA - END
273#endif
Note: See TracBrowser for help on using the repository browser.