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
|
/*@@
@file jacobi_wrapper.c
@date Tue Aug 24 12:50:07 1999
@author Gerd Lanfermann
@desc
The C wrapper, which calles the core Fortran routine, which
performs the actual solve.
We cannot derive the pointers to the GF data from the indeces in
Fortran. So we do this here in C and then pass the everything
over to the Fortran routine.
This wrapper is registers with the Elliptic solver registry
(not the Fortran file) , as coded up in ./CactusElliptic/EllBase
@enddesc
@@*/
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_FortranString.h"
void FORTRAN_NAME(sor_confmetric_core3d)
(_CCTK_C2F_PROTO(GH),
int *,
CCTK_REAL *,
int *,
CCTK_REAL *,
CCTK_REAL *,
CCTK_REAL *,
CCTK_REAL *,
CCTK_REAL *,
CCTK_REAL *,
CCTK_REAL *,
CCTK_REAL *,
CCTK_REAL *,
int *,
CCTK_REAL *,
CCTK_REAL *);
void FORTRAN_NAME(sor_flat_core3d)
(int *ierr,
_CCTK_C2F_PROTO(GH),
int *,
CCTK_REAL *,
int *,
CCTK_REAL *,
CCTK_REAL *,
int *,
CCTK_REAL *,
CCTK_REAL *);
/* We pass in the arguments that are neccessary for this class of elliptic eq.
this solver is intended to solve. See ./CactusElliptic/EllBase/src/ for the
classes of elliptic eq. */
void sor_confmetric(cGH *GH,
int *MetricPsiI,
int FieldIndex,
int MIndex,
int NIndex,
CCTK_REAL *AbsTol,
CCTK_REAL *RelTol)
{
CCTK_REAL *gxx=NULL, *gxy=NULL, *gxz=NULL;
CCTK_REAL *gyy=NULL, *gyz=NULL, *gzz=NULL;
CCTK_REAL *psi=NULL;
CCTK_REAL *Mlinear=NULL, *Nsources=NULL;
CCTK_REAL *Field =NULL;
int i;
int Mlinear_lsh[3], Nsource_lsh[3];
/* derive the metric data pointer from the index array. Note the ordering.
Also get datapointers to the field to solve for.
All of these are mandatory */
gxx = (CCTK_REAL*) CCTK_VarDataPtrI(GH, 0, MetricPsiI[0]);
gxy = (CCTK_REAL*) CCTK_VarDataPtrI(GH, 0, MetricPsiI[1]);
gxz = (CCTK_REAL*) CCTK_VarDataPtrI(GH, 0, MetricPsiI[2]);
gyy = (CCTK_REAL*) CCTK_VarDataPtrI(GH, 0, MetricPsiI[3]);
gyz = (CCTK_REAL*) CCTK_VarDataPtrI(GH, 0, MetricPsiI[4]);
gzz = (CCTK_REAL*) CCTK_VarDataPtrI(GH, 0, MetricPsiI[5]);
psi = (CCTK_REAL*) CCTK_VarDataPtrI(GH, 0, MetricPsiI[6]);
Field = (CCTK_REAL*) CCTK_VarDataPtrI(GH,0,FieldIndex);
if ((!gxx)||(!gxy)||(!gxz)||(!gyy)||(!gyz)||(!gzz)||(!psi)||(!Field))
{
CCTK_WARN(0,"SOR_WRAPPER: One of the metric data fields, or the GF to solve could not be found!");
}
/* derive the data pointer for the fields. the M/N fields are not
allocated (better: are of size 1), if the passed index is negative,
or we get back an empty GF of size 1 */
if (MIndex>=0) Mlinear = (CCTK_REAL*) CCTK_VarDataPtrI(GH,0,MIndex);
if (NIndex>=0) Nsources = (CCTK_REAL*) CCTK_VarDataPtrI(GH,0,NIndex);
/* we pass the size of M/N through to fortran, so F can
tell the difference between an allocated GF (Mlinear_lsh.eq.cctk_lsh) or
unallocated GF (Mlinear_lsh=1). maximal dimension is three. */
if (GH->cctk_dim>3)
CCTK_WARN(0,"This elliptic solver implementation does not do dimension>3!");
for (i=0;i<GH->cctk_dim;i++)
{
if((MIndex<0)) Mlinear_lsh[i]=1;
else Mlinear_lsh[i]=GH->cctk_lsh[i];
if((NIndex<0)) Nsource_lsh[i]=1;
else Nsource_lsh[i]=GH->cctk_lsh[i];
}
/* call the fortran routine */
FORTRAN_NAME(sor_confmetric_core3d)(_PASS_CCTK_C2F(GH),
Mlinear_lsh, Mlinear,
Nsource_lsh, Nsources,
gxx,gxy,gxz,gyy,gyz,gzz,psi,
Field, &FieldIndex, AbsTol, RelTol);
}
int sor_flat(cGH *GH,
int FieldIndex,
int MIndex,
int NIndex,
CCTK_REAL *AbsTol,
CCTK_REAL *RelTol)
{
int ierr;
CCTK_REAL *Mlinear=NULL, *Nsources=NULL;
CCTK_REAL *Field=NULL;
int i;
int Mlinear_lsh[3], Nsource_lsh[3];
Field = (CCTK_REAL*) CCTK_VarDataPtrI(GH,0,FieldIndex);
if (MIndex>0) Mlinear = (CCTK_REAL*) CCTK_VarDataPtrI(GH,0,MIndex);
if (NIndex>0) Nsources = (CCTK_REAL*) CCTK_VarDataPtrI(GH,0,NIndex);
if (GH->cctk_dim>3)
{
CCTK_WARN(0,"This elliptic solver implementation does not do dimension>3!");
}
for (i=0;i<GH->cctk_dim;i++)
{
if((MIndex<0)) Mlinear_lsh[i]=1;
else Mlinear_lsh[i]=GH->cctk_lsh[i];
if((NIndex<0)) Nsource_lsh[i]=1;
else
{
Nsource_lsh[i]=GH->cctk_lsh[i];
}
}
/* call the fortran routine */
FORTRAN_NAME(sor_flat_core3d)
(&ierr,
_PASS_CCTK_C2F(GH),
Mlinear_lsh,
Mlinear,
Nsource_lsh,
Nsources,
Field,
&FieldIndex,
AbsTol,
RelTol);
return ierr;
}
|