ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/SHAPES/SphereHarm.cpp
Revision: 1298
Committed: Thu Jun 24 15:47:52 2004 UTC (20 years ago) by chrisfen
File size: 6090 byte(s)
Log Message:
Shape file output format adjusted

File Contents

# Content
1 #include "SphereHarm.hpp"
2
3 SphereHarm::SphereHarm( int bandWidth ){
4 bw = bandWidth;
5
6 /*** ASSUMING WILL SEMINAIVE ALL ORDERS ***/
7 cutoff = bw;
8 size = 2*bw;
9
10 /* allocate memory */
11 rdata = (double *) fftw_malloc(sizeof(double) * (size * size));
12 idata = (double *) fftw_malloc(sizeof(double) * (size * size));
13 rcoeffs = (double *) fftw_malloc(sizeof(double) * (bw * bw));
14 icoeffs = (double *) fftw_malloc(sizeof(double) * (bw * bw));
15 weights = (double *) fftw_malloc(sizeof(double) * 4 * bw);
16 seminaive_naive_tablespace =
17 (double *) fftw_malloc(sizeof(double) *
18 (Reduced_Naive_TableSize(bw,cutoff) +
19 Reduced_SpharmonicTableSize(bw,cutoff)));
20 workspace = (double *) fftw_malloc(sizeof(double) *
21 ((8 * (bw*bw)) +
22 (7 * bw)));
23
24
25 /****
26 At this point, check to see if all the memory has been
27 allocated. If it has not, there's no point in going further.
28 ****/
29
30 if ( (rdata == NULL) || (idata == NULL) ||
31 (rcoeffs == NULL) || (icoeffs == NULL) ||
32 (seminaive_naive_tablespace == NULL) ||
33 (workspace == NULL) )
34 {
35 perror("Error in allocating memory");
36 exit( 1 ) ;
37 }
38
39 //precompute the Legendres
40 fprintf(stdout,"Precomputing the Legendres...\n");
41 seminaive_naive_table = SemiNaive_Naive_Pml_Table( bw, cutoff,
42 seminaive_naive_tablespace,
43 workspace );
44
45 //construct fftw plans using the GURU interface
46 /* forward DCT */
47 dctPlan = fftw_plan_r2r_1d( 2*bw, weights, rdata,
48 FFTW_REDFT10, FFTW_ESTIMATE ) ;
49
50 /*
51 fftw "preamble" ;
52 note that this plan places the output in a transposed array
53 */
54 rank = 1 ;
55 dims[0].n = 2*bw ;
56 dims[0].is = 1 ;
57 dims[0].os = 2*bw ;
58 howmany_rank = 1 ;
59 howmany_dims[0].n = 2*bw ;
60 howmany_dims[0].is = 2*bw ;
61 howmany_dims[0].os = 1 ;
62
63 /* forward fft */
64 fftPlan = fftw_plan_guru_split_dft( rank, dims,
65 howmany_rank, howmany_dims,
66 rdata, idata,
67 workspace, workspace+(4*bw*bw),
68 FFTW_ESTIMATE );
69
70 //make the weights
71 makeweights( bw, weights );
72 }
73
74 SphereHarm::~SphereHarm(){
75 //free up memory
76 fftw_destroy_plan( fftPlan );
77 fftw_destroy_plan( dctPlan );
78
79 fftw_free(workspace);
80 fftw_free(seminaive_naive_table);
81 fftw_free(seminaive_naive_tablespace);
82 fftw_free(weights);
83 fftw_free(icoeffs);
84 fftw_free(rcoeffs);
85 fftw_free(idata);
86 fftw_free(rdata);
87 }
88
89 void SphereHarm::doTransforms(vector<double> gridData){
90 int i;
91
92 //load the data
93 for (i=0; i<size*size; i++){
94 rdata[i] = gridData[i];
95 //our data is all real, so load the imaginary part with zeros
96 idata[i] = 0.0;
97 }
98
99 //do the forward spherical transform
100 FST_semi_memo(rdata, idata,
101 rcoeffs, icoeffs,
102 bw,
103 seminaive_naive_table,
104 workspace,
105 0,
106 cutoff,
107 &dctPlan,
108 &fftPlan,
109 weights );
110 }
111
112 void SphereHarm::printShapesFileStart(char name[200], char particle[80],
113 double mass, double momInert[3][3]){
114 ofstream shapes(name);
115 shapes << "begin ShapeInfo\n";
116 shapes << "#name\t\tmass\tI_xx\tI_yy\tI_zz\n";
117 shapes << particle << "\t" << mass << "\t" << momInert[0][0] << "\t"
118 << momInert[1][1] << "\t" << momInert[2][2] << "\n";
119 shapes << "end ShapeInfo\n";
120 }
121
122 void SphereHarm::printToShapesFile(char name[200], int index){
123 ofstream shapes(name, ios::app);
124
125 biggest = 0.0;
126 nfuncs = 0;
127 for ( l = 0 ; l < bw ; l++ ) {
128 for (m = 0; m < l+1; m++) {
129 dummy1 = seanindex(m, l, bw);
130 dummy2 = seanindex(-m, l, bw);
131
132 if (m == 0) {
133 cm = rcoeffs[dummy1];
134 sm = icoeffs[dummy1];
135 } else {
136 cm = pow(-1.0,(double)m)*rcoeffs[dummy1] + rcoeffs[dummy2];
137 sm = pow(-1.0,(double)m)*icoeffs[dummy1] - icoeffs[dummy2];
138 }
139
140 if (fabs(cm) > biggest) biggest = fabs(cm);
141 if (fabs(sm) > biggest) biggest = fabs(sm);
142 }
143 }
144 for ( l = 0 ; l < bw ; l++ ) {
145 for (m = 0; m < l+1; m++) {
146 dummy1 = seanindex(m, l, bw);
147 dummy2 = seanindex(-m, l, bw);
148
149 if (m == 0) {
150 cm = rcoeffs[dummy1];
151 sm = icoeffs[dummy1];
152 } else {
153 cm = pow(-1.0,(double)m)*rcoeffs[dummy1] + rcoeffs[dummy2];
154 sm = pow(-1.0,(double)m)*icoeffs[dummy1] - icoeffs[dummy2];
155 }
156
157 if (fabs(cm) > 0.01 * biggest) nfuncs++;
158 if (fabs(sm) > 0.01 * biggest) nfuncs++;
159 }
160 }
161
162 switch(index){
163 case 0:{
164 shapes << "\nbegin ContactFunctions\n";
165 shapes << "#l\tm\tsin or cos\tcoeff (Ang)\n";
166 }; break;
167 case 1:{
168 shapes << "\nbegin RangeFunctions\n";
169 shapes << "#l\tm\tsin or cos\tcoeff (Ang)\n";
170 }; break;
171 case 2:{
172 shapes << "\nbegin StrengthFunctions\n";
173 shapes << "#l\tm\tsin or cos\tcoeff (kcal/mol)\n";
174 }; break;
175 }
176
177 for ( l = 0 ; l < bw ; l++ ) {
178 for (m = 0; m < l+1; m++) {
179 dummy1 = seanindex(m, l, bw);
180 dummy2 = seanindex(-m, l, bw);
181
182 if (m == 0) {
183 cm = rcoeffs[dummy1];
184 sm = icoeffs[dummy1];
185 } else {
186 cm = pow(-1.0,(double)m)*rcoeffs[dummy1] + rcoeffs[dummy2];
187 sm = pow(-1.0,(double)m)*icoeffs[dummy1] - icoeffs[dummy2];
188 }
189
190 if (fabs(cm) > 0.01 * biggest)
191 shapes << l << "\t" << m << "\tcos\t\t" << cm << "\n";
192 if (fabs(sm) > 0.01 * biggest)
193 shapes << l << "\t" << m << "\tsin\t\t" << sm << "\n";
194 }
195 }
196 switch(index){
197 case 0:{
198 shapes << "end ContactFunctions\n";
199 }; break;
200 case 1:{
201 shapes << "end RangeFunctions\n";
202 }; break;
203 case 2:{
204 shapes << "end StrengthFunctions\n";
205 }; break;
206 }
207 }
208