aboutsummaryrefslogtreecommitdiff
path: root/src/elliptic/Jacobian.cc
blob: 58009ec1c2d6f4da689073e5478ee711fff9eb9d (plain)
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
// Jacobian.cc -- data structures for the Jacobian matrix
// $Id$

//
// Jacobian::Jacobian
//
// dense_Jacobian::dense_Jacobian
// dense_Jacobian::~dense_Jacobian
// dense_Jacobian::zero
// dense_Jacobian::zero_row
// dense_Jacobian::solve_linear_system
//

#include <stdio.h>
using std::fopen;
using std::printf;
using std::fprintf;
using std::fclose;
using std::FILE;
#include <assert.h>
#include <math.h>
#include <vector>

#include "util_Table.h"
#include "cctk.h"
#include "cctk_Arguments.h"

#include "stdc.h"
#include "config.hh"
#include "../jtutil/util.hh"
#include "../jtutil/array.hh"
#include "../jtutil/cpm_map.hh"
#include "../jtutil/linear_map.hh"
using jtutil::error_exit;

#include "../util/coords.hh"
#include "../util/grid.hh"
#include "../util/fd_grid.hh"
#include "../util/patch.hh"
#include "../util/patch_edge.hh"
#include "../util/patch_interp.hh"
#include "../util/ghost_zone.hh"
#include "../util/patch_system.hh"

#include "Jacobian.hh"
//#include "lapack.h"
//**************************************
/* lapack.h -- C/C++ prototypes for (some) LAPACK routines */
/* $Id$ */

/*
 * prerequisites:
 *	"cctk.h"
 *	"config.hh"	// for "integer" = Fortran integer
 */

#ifdef __cplusplus
extern "C" {
#endif

void CCTK_FCALL
  CCTK_FNAME(sgesv)(const integer* N, const integer* NRHS,
		    float A[], const int* LDA,
		    integer IPIV[],
		    float B[], const integer* LDB, integer* info);
void CCTK_FCALL
  CCTK_FNAME(dgesv)(const integer* N, const integer* NRHS,
		    double A[], const int* LDA,
		    integer IPIV[],
		    double B[], const integer* LDB, integer* info);

#ifdef __cplusplus
           };	/* extern "C" */
#endif
//**************************************

//******************************************************************************
//******************************************************************************
//******************************************************************************

//
// This function constructs a  Jacobian  object.
//
Jacobian::Jacobian(patch_system& ps)
	: ps_(ps),
	  NN_(ps.N_grid_points())
{ }

//******************************************************************************
//******************************************************************************
//******************************************************************************

//
// This function constructs a  dense_Jacobian  object.
//
dense_Jacobian::dense_Jacobian(patch_system& ps)
	: Jacobian(ps),
	  matrix_(0,NN()-1, 0,NN()-1),
	  pivot_(new integer[NN()])
{ }

//******************************************************************************

//
// THis function destroys a  dense_Jacobian  object.
//
dense_Jacobian::~dense_Jacobian()
{
delete[] pivot_;
}

//******************************************************************************

//
// This function zeros a  dense_Jacobian  object.
//
void dense_Jacobian::zero()
{
	for (int JJ = 0 ; JJ < NN() ; ++JJ)
	{
		for (int II = 0 ; II < NN() ; ++II)
		{
		matrix_(JJ,II) = 0.0;
		}
	}
}

//******************************************************************************

//
// This function zeros a single row of a  dense_Jacobian  object.
//
void dense_Jacobian::zero_row(int II)
{
	for (int JJ = 0 ; JJ < NN() ; ++JJ)
	{
	matrix_(JJ,II) = 0.0;
	}
}

//******************************************************************************

//
// This function solves the linear system J.x = rhs, with rhs and x
// being nominal-grid gridfns.  The computation is done using the LAPACK
// [sd]gesv() subroutines; which modify the Jacobian matrix J in-place
// for the LU decomposition.
//
void dense_Jacobian::solve_linear_system(int rhs_gfn, int x_gfn)
{
const fp *rhs = my_patch_system().gridfn_data(rhs_gfn);
fp       *x   = my_patch_system().gridfn_data(x_gfn);

//
// [sd]gesv() use an "in out" design, where the same argument is used for
// both rhs and x ==> first copy rhs to x so we can pass that to [sd]gesv()
//
	for (int II = 0 ; II < NN() ; ++II)
	{
	x[II] = rhs[II];
	}

integer N = NN();
integer NRHS = 1;
integer info;

#ifdef FP_IS_FLOAT
  CCTK_FNAME(sgesv)(&N, &NRHS, matrix_.data_array(), &N, pivot_, x, &N, &info);
#endif
#ifdef FP_IS_DOUBLE
  CCTK_FNAME(dgesv)(&N, &NRHS, matrix_.data_array(), &N, pivot_, x, &N, &info);
#endif

if (info != 0)
   then error_exit(ERROR_EXIT,
"\n"
"***** dense_Jacobian::solve_linear_system(rhs_gfn=%d, x_gfn=%d):\n"
"        error return info=%d from [sd]gesv() LAPACK routine!\n"
,
		   rhs_gfn, x_gfn,
		   int(info));					/*NOTREACHED*/
}