Epetra Package Browser (Single Doxygen Collection)
Development
Toggle main menu visibility
Loading...
Searching...
No Matches
src
Epetra_SerialSymDenseMatrix.cpp
Go to the documentation of this file.
1
2
//@HEADER
3
// ************************************************************************
4
//
5
// Epetra: Linear Algebra Services Package
6
// Copyright 2011 Sandia Corporation
7
//
8
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9
// the U.S. Government retains certain rights in this software.
10
//
11
// Redistribution and use in source and binary forms, with or without
12
// modification, are permitted provided that the following conditions are
13
// met:
14
//
15
// 1. Redistributions of source code must retain the above copyright
16
// notice, this list of conditions and the following disclaimer.
17
//
18
// 2. Redistributions in binary form must reproduce the above copyright
19
// notice, this list of conditions and the following disclaimer in the
20
// documentation and/or other materials provided with the distribution.
21
//
22
// 3. Neither the name of the Corporation nor the names of the
23
// contributors may be used to endorse or promote products derived from
24
// this software without specific prior written permission.
25
//
26
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37
//
38
// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39
//
40
// ************************************************************************
41
//@HEADER
42
43
#include "
Epetra_SerialSymDenseMatrix.h
"
44
//=============================================================================
45
Epetra_SerialSymDenseMatrix::Epetra_SerialSymDenseMatrix
(
void
)
46
:
Epetra_SerialDenseMatrix
(),
47
Upper_
(false),
48
UPLO_
(
'L'
)
49
50
{
51
}
52
//=============================================================================
53
Epetra_SerialSymDenseMatrix::Epetra_SerialSymDenseMatrix
(
Epetra_DataAccess
CV_in,
double
*A_in,
int
LDA_in,
int
NumRowsCols)
54
:
Epetra_SerialDenseMatrix
(CV_in, A_in, LDA_in, NumRowsCols, NumRowsCols),
55
Upper_
(false),
56
UPLO_
(
'L'
)
57
58
{
59
}
60
//=============================================================================
61
Epetra_SerialSymDenseMatrix::Epetra_SerialSymDenseMatrix
(
const
Epetra_SerialSymDenseMatrix
& Source)
62
:
Epetra_SerialDenseMatrix
(Source),
63
Upper_
(Source.
Upper_
),
64
UPLO_
(Source.
UPLO_
)
65
{
66
}
67
//=============================================================================
68
Epetra_SerialSymDenseMatrix::~Epetra_SerialSymDenseMatrix
()
69
{
70
}
71
//=============================================================================
72
void
Epetra_SerialSymDenseMatrix::CopyUPLOMat
(
bool
Upper_in,
double
* A_in,
int
LDA_in,
int
NumRows) {
73
74
int
i, j;
75
double
* ptr1;
76
double
* ptr2;
77
78
if
(Upper_in) {
79
for
(j=1; j<NumRows; j++) {
80
ptr1 = A_in + j;
81
ptr2 = A_in + j*LDA_in;
82
for
(i=0; i<j; i++) {
83
*ptr1 = *ptr2++;
84
ptr1+=LDA_in;
85
}
86
}
87
}
88
else
{
89
for
(i=1; i<NumRows; i++) {
90
ptr1 = A_in + i;
91
ptr2 = A_in + i*LDA_in;
92
for
(j=0; j<i; j++) {
93
*ptr2++ = *ptr1;
94
ptr1+=LDA_in;
95
}
96
}
97
}
98
}
99
//=============================================================================
100
double
Epetra_SerialSymDenseMatrix::NormOne
(
void
)
const
{
101
102
return
(
Epetra_SerialSymDenseMatrix::NormInf
());
103
}
104
//=============================================================================
105
double
Epetra_SerialSymDenseMatrix::NormInf
(
void
)
const
{
106
107
int
i, j;
108
109
double
anorm = 0.0;
110
double
* ptr;
111
112
if
(!
Upper
()) {
113
for
(j=0; j<
N_
; j++) {
114
double
sum = 0.0;
115
ptr =
A_
+ j + j*
LDA_
;
116
for
(i=j; i<
N_
; i++) sum += std::abs(*ptr++);
117
ptr =
A_
+ j;
118
for
(i=0; i<j; i++) {
119
sum += std::abs(*ptr);
120
ptr +=
LDA_
;
121
}
122
anorm =
EPETRA_MAX
(anorm, sum);
123
}
124
}
125
else
{
126
for
(j=0; j<
N_
; j++) {
127
double
sum = 0.0;
128
ptr =
A_
+ j*
LDA_
;
129
for
(i=0; i<j; i++) sum += std::abs(*ptr++);
130
ptr =
A_
+ j + j*
LDA_
;
131
for
(i=j; i<
N_
; i++) {
132
sum += std::abs(*ptr);
133
ptr +=
LDA_
;
134
}
135
anorm =
EPETRA_MAX
(anorm, sum);
136
}
137
}
138
UpdateFlops
(
N_
*
N_
);
139
return
(anorm);
140
}
141
142
//=============================================================================
143
int
Epetra_SerialSymDenseMatrix::Scale
(
double
ScalarA) {
144
145
int
i, j;
146
147
double
* ptr;
148
149
if
(!
Upper
()) {
150
for
(j=0; j<
N_
; j++) {
151
ptr =
A_
+ j + j*
LDA_
;
152
for
(i=j; i<
N_
; i++) {*ptr = *ptr * ScalarA; ptr++;}
153
}
154
}
155
else
{
156
for
(j=0; j<
N_
; j++) {
157
ptr =
A_
+ j*
LDA_
;
158
for
(i=0; i<j; i++) {*ptr = *ptr * ScalarA; ptr++;}
159
}
160
}
161
UpdateFlops
(
N_
*(
N_
+1)/2);
162
return
(0);
163
}
164
EPETRA_MAX
#define EPETRA_MAX(x, y)
Definition
Epetra_ConfigDefs.h:62
Epetra_DataAccess
Epetra_DataAccess
Definition
Epetra_DataAccess.h:55
Epetra_SerialSymDenseMatrix.h
Epetra_CompObject::UpdateFlops
void UpdateFlops(int Flops_in) const
Increment Flop count for this object.
Definition
Epetra_CompObject.h:99
Epetra_SerialDenseMatrix::LDA_
int LDA_
Definition
Epetra_SerialDenseMatrix.h:491
Epetra_SerialDenseMatrix::A_
double * A_
Definition
Epetra_SerialDenseMatrix.h:492
Epetra_SerialDenseMatrix::Epetra_SerialDenseMatrix
Epetra_SerialDenseMatrix(bool set_object_label=true)
Default constructor; defines a zero size object.
Definition
Epetra_SerialDenseMatrix.cpp:49
Epetra_SerialDenseMatrix::N_
int N_
Definition
Epetra_SerialDenseMatrix.h:480
Epetra_SerialSymDenseMatrix::UPLO_
char UPLO_
Definition
Epetra_SerialSymDenseMatrix.h:277
Epetra_SerialSymDenseMatrix::CopyUPLOMat
void CopyUPLOMat(bool Upper, double *A, int LDA, int NumRows)
Definition
Epetra_SerialSymDenseMatrix.cpp:72
Epetra_SerialSymDenseMatrix::Epetra_SerialSymDenseMatrix
Epetra_SerialSymDenseMatrix(void)
Default constructor; defines a zero size object.
Definition
Epetra_SerialSymDenseMatrix.cpp:45
Epetra_SerialSymDenseMatrix::~Epetra_SerialSymDenseMatrix
virtual ~Epetra_SerialSymDenseMatrix()
Epetra_SerialSymDenseMatrix destructor.
Definition
Epetra_SerialSymDenseMatrix.cpp:68
Epetra_SerialSymDenseMatrix::Upper
bool Upper() const
Returns true if upper triangle of this matrix has and will be used.
Definition
Epetra_SerialSymDenseMatrix.h:225
Epetra_SerialSymDenseMatrix::NormInf
double NormInf() const
Computes the Infinity-Norm of the this matrix.
Definition
Epetra_SerialSymDenseMatrix.cpp:105
Epetra_SerialSymDenseMatrix::Upper_
bool Upper_
Definition
Epetra_SerialSymDenseMatrix.h:275
Epetra_SerialSymDenseMatrix::NormOne
double NormOne() const
Computes the 1-Norm of the this matrix.
Definition
Epetra_SerialSymDenseMatrix.cpp:100
Epetra_SerialSymDenseMatrix::Scale
int Scale(double ScalarA)
Inplace scalar-matrix product A = a A.
Definition
Epetra_SerialSymDenseMatrix.cpp:143
Generated by
1.17.0