e3bf30522c36ea94d6c348766762da27092722d1
[alexxy/gromacs.git] / src / gmxlib / gmx_blas / dsyr2k.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012, by the GROMACS development team, led by
5  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
6  * others, as listed in the AUTHORS file in the top-level source
7  * directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 #include <math.h>
36 #include <ctype.h>
37
38
39 #include <types/simple.h>
40 #include "gmx_blas.h"
41
42 void
43 F77_FUNC(dsyr2k,DSYR2K)(const char *uplo, 
44         const char *trans,
45         int *n__,
46         int *k__,
47         double *alpha__,
48         double *a,
49         int *lda__,
50         double *b,
51         int *ldb__,
52         double *beta__,
53         double *c,
54         int *ldc__)
55 {
56   char ch1,ch2;
57   int nrowa;
58   int i,j,l;
59   double temp1,temp2;
60
61   
62   int n = *n__;
63   int k = *k__;
64   int lda = *lda__;
65   int ldb = *ldb__;
66   int ldc = *ldc__;
67   
68   double alpha = *alpha__;
69   double beta  = *beta__;
70   
71   ch1 = toupper(*uplo);
72   ch2 = toupper(*trans);
73
74   if(ch2 == 'N')
75     nrowa = n;
76   else
77     nrowa = k;
78
79   if(n==0 || ( ( fabs(alpha)<GMX_DOUBLE_MIN || k==0 ) && fabs(beta-1.0)<GMX_DOUBLE_EPS))
80     return;
81
82   if(fabs(alpha)<GMX_DOUBLE_MIN ) {
83     if(ch1=='U') {
84       if(fabs(beta)<GMX_DOUBLE_MIN) 
85         for(j=1;j<=n;j++) 
86           for(i=1;i<=j;i++)
87             c[(j-1)*(ldc)+(i-1)] = 0.0;
88       else
89         for(j=1;j<=n;j++) 
90           for(i=1;i<=j;i++)
91             c[(j-1)*(ldc)+(i-1)] *= beta;
92     } else {
93       /* lower */
94       if(fabs(beta)<GMX_DOUBLE_MIN) 
95         for(j=1;j<=n;j++) 
96           for(i=j;i<=n;i++)
97             c[(j-1)*(ldc)+(i-1)] = 0.0;
98       else
99         for(j=1;j<=n;j++) 
100           for(i=j;i<=n;i++)
101             c[(j-1)*(ldc)+(i-1)] *= beta;
102     }
103     return;
104   }
105
106   if(ch2=='N') {
107     if(ch1=='U') {
108       for(j=1;j<=n;j++) {
109         if(fabs(beta)<GMX_DOUBLE_MIN)
110           for(i=1;i<=j;i++)
111              c[(j-1)*(ldc)+(i-1)] = 0.0;
112         else if(fabs(beta-1.0)>GMX_DOUBLE_EPS)
113           for(i=1;i<=j;i++)
114             c[(j-1)*(ldc)+(i-1)] *= beta;
115         for(l=1;l<=k;l++) {
116           if( fabs(a[(l-1)*(lda)+(j-1)])>GMX_DOUBLE_MIN ||
117               fabs(b[(l-1)*(ldb)+(j-1)])>GMX_DOUBLE_MIN) {
118             temp1 = alpha * b[(l-1)*(ldb)+(j-1)];
119             temp2 = alpha * a[(l-1)*(lda)+(j-1)];
120             for(i=1;i<=j;i++)
121               c[(j-1)*(ldc)+(i-1)] += 
122                 a[(l-1)*(lda)+(i-1)] * temp1 + 
123                 b[(l-1)*(ldb)+(i-1)] * temp2;
124           }
125         }
126       }
127     } else {
128       /* lower */
129       for(j=1;j<=n;j++) {
130         if(fabs(beta)<GMX_DOUBLE_MIN)
131           for(i=j;i<=n;i++)
132             c[(j-1)*(ldc)+(i-1)] = 0.0;
133         else if(fabs(beta-1.0)>GMX_DOUBLE_EPS)
134           for(i=j;i<=n;i++)
135             c[(j-1)*(ldc)+(i-1)] *= beta;
136         for(l=1;l<=k;l++) {
137           if( fabs(a[(l-1)*(lda)+(j-1)])>GMX_DOUBLE_MIN ||
138               fabs(b[(l-1)*(ldb)+(j-1)])>GMX_DOUBLE_MIN) {
139             temp1 = alpha * b[(l-1)*(ldb)+(j-1)];
140             temp2 = alpha * a[(l-1)*(lda)+(j-1)];
141             for(i=j;i<=n;i++)
142               c[(j-1)*(ldc)+(i-1)] += 
143                 a[(l-1)*(lda)+(i-1)] * temp1 + 
144                 b[(l-1)*(ldb)+(i-1)] * temp2;
145           }
146         }
147       }
148     }
149   } else {
150     /* transpose */
151     if(ch1=='U') {
152       for(j=1;j<=n;j++) 
153         for(i=1;i<=j;i++) {
154           temp1 = 0.0;
155           temp2 = 0.0;
156           for (l=1;l<=k;l++) {
157              temp1 += a[(i-1)*(lda)+(l-1)] * b[(j-1)*(ldb)+(l-1)];
158              temp2 += b[(i-1)*(ldb)+(l-1)] * a[(j-1)*(lda)+(l-1)];
159           }
160           if(fabs(beta)<GMX_DOUBLE_MIN)
161             c[(j-1)*(ldc)+(i-1)] = alpha * (temp1 + temp2);
162           else
163             c[(j-1)*(ldc)+(i-1)] = beta * c[(j-1)*(ldc)+(i-1)] +
164               alpha * (temp1 + temp2);
165         }
166     } else {
167       /* lower */
168       for(j=1;j<=n;j++) 
169         for(i=j;i<=n;i++) {
170           temp1 = 0.0;
171           temp2 = 0.0;
172           for (l=1;l<=k;l++) {
173              temp1 += a[(i-1)*(lda)+(l-1)] * b[(j-1)*(ldb)+(l-1)];
174              temp2 += b[(i-1)*(ldb)+(l-1)] * a[(j-1)*(lda)+(l-1)];
175           }
176           if(fabs(beta)<GMX_DOUBLE_MIN)
177             c[(j-1)*(ldc)+(i-1)] = alpha * (temp1 + temp2);
178           else
179             c[(j-1)*(ldc)+(i-1)] = beta * c[(j-1)*(ldc)+(i-1)] +
180               alpha * (temp1 + temp2);
181         }
182     }
183   }
184   return;
185 }