aboutsummaryrefslogtreecommitdiff
path: root/src/D3_to_D2.F
blob: a7d22b2754ba79b79d65a7681e437098ae5efdb9 (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
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
#include "cctk.h"

c     ========================================================================

      SUBROUTINE D3_to_D2(cctkGH,conformal_state,do_ADMmass,do_momentum,do_spin,
     &     Psi_power,origin,myproc,interpolation_order,Dx,Dy,Dz,Psi,
     &     g00,gxx,gxy,gxz,gyy,gyz,gzz,hxx,hxy,hxz,hyy,hyz,hzz,
     &     x,y,z,eta,Nt,Np,theta,phi,
     &     Psis,g00s,gxxs,gxys,gxzs,gyys,gyzs,gzzs,dPsis,dgxxs,dgxys,dgxzs,
     &     dgyys,dgyzs,dgzzs,ADMmass_int1,ADMmass_int2,
     &     momentum_int1,momentum_int2,momentum_int3,
     &     spin_int1,spin_int2,spin_int3,Extract_temp3d)
 
c     ------------------------------------------------------------------------
c
c     Project the 3-metric and its 1st radial derivatives onto a 
c     2-sphere.
c
c     ------------------------------------------------------------------------

      USE met_rad_der_int
      USE ADMmass_integrand3D_int
      USE momentum_integrand3D_int
      USE spin_integrand3D_int

      IMPLICIT NONE

c     Input variables

      CCTK_POINTER :: cctkGH

      INTEGER,INTENT(IN) :: 
     &     conformal_state,myproc,Psi_power
      CCTK_INT, INTENT(IN) ::
     &     Nt,Np,do_momentum,do_spin,interpolation_order
      CCTK_REAL,INTENT(IN) :: 
     &     origin(3),Dx,Dy,Dz,eta
      CCTK_REAL,INTENT(IN),DIMENSION(:) :: 
     &     theta,phi,x,y,z
      CCTK_REAL,INTENT(IN),DIMENSION(:,:,:) ::  
     &     Psi,g00,gxx,gxy,gxz,gyy,gyz,gzz,
     &     hxx,hxy,hxz,hyy,hyz,hzz
      CCTK_REAL,INTENT(INOUT),DIMENSION(:,:,:) ::
     &     Extract_temp3d
      LOGICAL ::
     &     do_ADMmass(2)

c     Output variables

      CCTK_REAL,INTENT(OUT),DIMENSION(:,:) ::     
     &     Psis,g00s,gxxs,gxys,gxzs,gyys,
     &     gyzs,gzzs,dPsis,dgxxs,dgxys,dgxzs,dgyys,dgyzs,dgzzs,
     &     ADMmass_int1,ADMmass_int2,
     &     momentum_int1,momentum_int2,momentum_int3,
     &     spin_int1,spin_int2,spin_int3


c     Local variables, passed on

      LOGICAL :: 
     &    err_flag
      INTEGER :: 
     &    iorder,npoints,Nx,Ny,Nz
      INTEGER,DIMENSION(Nt,Np) :: 
     &    ib,jb,kb
      CCTK_REAL,DIMENSION(Nt,Np) :: 
     &    xs,ys,zs,ux,uy,uz,xb,yb,zb

c     Local variables, here only

      INTEGER :: 
     &    i,j,interp_handle,coord_system_handle,ierror
      INTEGER,   DIMENSION(8)   :: in_array_indices


c     ------------------------------------------------------------------------
      
c     Initial Stuff
c     -------------
      Nx = SIZE(x)
      Ny = SIZE(y)
      Nz = SIZE(z)
      

c     Compute Cartesian coordinates on the surface
c     --------------------------------------------
      DO j = 1, Np
         DO i = 1, Nt

            xs(i,j) = origin(1)+eta*SIN(theta(i))*COS(phi(j))
            ys(i,j) = origin(2)+eta*SIN(theta(i))*SIN(phi(j))
            zs(i,j) = origin(3)+eta*COS(theta(i))

         ENDDO
      ENDDO


c     Only do interpolation on one processor
c     --------------------------------------
      SELECT CASE (myproc) 

      CASE (0)
         npoints = Np*Nt

      CASE DEFAULT
         npoints = 0

      END SELECT 


c     Get the interpolator and coordinate system handle
c     -------------------------------------------------
      interp_handle = -1
      coord_system_handle = -1

      if (interpolation_order .eq. 1) then
        call CCTK_InterpHandle (interp_handle, "first-order uniform cartesian")
      else if (interpolation_order .eq. 2) then
        call CCTK_InterpHandle (interp_handle, "second-order uniform cartesian")
      else if (interpolation_order .eq. 3) then
        call CCTK_InterpHandle (interp_handle, "third-order uniform cartesian")
      endif

      call CCTK_CoordSystemHandle (coord_system_handle, "cart3d")

      if (interp_handle .lt. 0 .or. coord_system_handle .lt. 0) then
        call CCTK_WARN (0, "Couldn't get handles for interpolation operator and/or coordinate system")
      endif

c     Get indices of GFs to interpolate
c     ---------------------------------
      call CCTK_VarIndex(in_array_indices(1), "staticconformal::psi")
      call CCTK_VarIndex(in_array_indices(2), "extract::g00")
      call CCTK_VarIndex(in_array_indices(3), "admbase::gxx")
      call CCTK_VarIndex(in_array_indices(4), "admbase::gxy")
      call CCTK_VarIndex(in_array_indices(5), "admbase::gxz")
      call CCTK_VarIndex(in_array_indices(6), "admbase::gyy")
      call CCTK_VarIndex(in_array_indices(7), "admbase::gyz")
      call CCTK_VarIndex(in_array_indices(8), "admbase::gzz")


c     Project un-physical metric and conformal factor  onto sphere
c     ------------------------------------------------------------
      if (conformal_state > 0) then
        call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                    npoints, 8, 8,
     $                    xs, ys, zs,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL,
     $                    in_array_indices(1), in_array_indices(2),
     $                    in_array_indices(3), in_array_indices(4),
     $                    in_array_indices(5), in_array_indices(6),
     $                    in_array_indices(7), in_array_indices(8),
     $                    Psis, g00s, gxxs, gxys, gxzs, gyys, gyzs, gzzs,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL)
      else
        call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                    npoints, 7, 7,
     $                    xs, ys, zs,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL,
     $                    in_array_indices(2),
     $                    in_array_indices(3), in_array_indices(4),
     $                    in_array_indices(5), in_array_indices(6),
     $                    in_array_indices(7), in_array_indices(8),
     $                    g00s, gxxs, gxys, gxzs, gyys, gyzs, gzzs,
     $                    CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL)
      end if


c     Calculate radial derivatives and project onto sphere
c     ----------------------------------------------------
      call CCTK_VarIndex(in_array_indices(1), "extract::temp3d")

      if (conformal_state > 0) then
         CALL met_rad_der(origin,Dx,Dy,Dz,x,y,z,Psi,Extract_temp3d)
         CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")   
         call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $        npoints, 1, 1,
     $        xs, ys, zs,
     $        CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $        CCTK_VARIABLE_REAL,
     $        in_array_indices(1),
     $        dPsis,
     $        CCTK_VARIABLE_REAL)
      end if

      CALL met_rad_der(origin,Dx,Dy,Dz,x,y,z,gxx,Extract_temp3d)
      CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")
      call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                    npoints, 1, 1,
     $                    xs, ys, zs,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL,
     $                    in_array_indices(1),
     $                    dgxxs,
     $                    CCTK_VARIABLE_REAL)

      CALL met_rad_der(origin,Dx,Dy,Dz,x,y,z,gxy,Extract_temp3d)
      CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")
      call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                    npoints, 1, 1,
     $                    xs, ys, zs,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL,
     $                    in_array_indices(1),
     $                    dgxys,
     $                    CCTK_VARIABLE_REAL)

      CALL met_rad_der(origin,Dx,Dy,Dz,x,y,z,gxz,Extract_temp3d)
      CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")
      call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                    npoints, 1, 1,
     $                    xs, ys, zs,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL,
     $                    in_array_indices(1),
     $                    dgxzs,
     $                    CCTK_VARIABLE_REAL)

      CALL met_rad_der(origin,Dx,Dy,Dz,x,y,z,gyy,Extract_temp3d)
      CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")
      call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                    npoints, 1, 1,
     $                    xs, ys, zs,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL,
     $                    in_array_indices(1),
     $                    dgyys,
     $                    CCTK_VARIABLE_REAL)

      CALL met_rad_der(origin,Dx,Dy,Dz,x,y,z,gyz,Extract_temp3d)
      CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")
      call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                    npoints, 1, 1,
     $                    xs, ys, zs,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL,
     $                    in_array_indices(1),
     $                    dgyzs,
     $                    CCTK_VARIABLE_REAL)

      CALL met_rad_der(origin,Dx,Dy,Dz,x,y,z,gzz,Extract_temp3d)
      CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")
      call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                    npoints, 1, 1,
     $                    xs, ys, zs,
     $                    CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                    CCTK_VARIABLE_REAL,
     $                    in_array_indices(1),
     $                    dgzzs,
     $                    CCTK_VARIABLE_REAL)


c     Calculate integrands for ADM masses
c     -----------------------------------
c     Standard equation
      IF (do_ADMmass(1)) THEN
         CALL ADMmass_integrand3D(origin,Dx,Dy,Dz,x,y,z,gxx,gxy,
     &        gxz,gyy,gyz,gzz,Extract_temp3d,Psi,Psi_power,conformal_state)
         CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")
         call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                       npoints, 1, 1,
     $                       xs, ys, zs,
     $                       CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                       CCTK_VARIABLE_REAL,
     $                       in_array_indices(1),
     $                       ADMmass_int1,
     $                       CCTK_VARIABLE_REAL)

      END IF

c     Conformal equation
      IF (do_ADMmass(2)) THEN
             ADMmass_int2 = -eta**2/2D0/3.1416D0*dPsis
      ENDIF	

c     Calculate integrands for momentum
c     ---------------------------------
      IF (do_momentum==1) THEN

         CALL momentum_integrand3D(origin,Dx,Dy,Dz,x,y,z,
     &        gxx,gxy,gxz,gyy,gyz,gzz,hxx,hxy,hxz,hyy,hyz,hzz,
     &        Extract_temp3d,Psi,Psi_power,conformal_state)
         CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")
         call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                       npoints, 1, 1,
     $                       xs, ys, zs,
     $                       CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                       CCTK_VARIABLE_REAL,
     $                       in_array_indices(1),
     $                       momentum_int1,
     $                       CCTK_VARIABLE_REAL)

         CALL momentum_integrand3D(origin,Dx,Dy,Dz,x,y,z,
     &        gxx,gxy,gxz,gyy,gyz,gzz,hxx,hxy,hxz,hyy,hyz,hzz,
     &        Extract_temp3d,Psi,Psi_power,conformal_state)
         CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")
         call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                       npoints, 1, 1,
     $                       xs, ys, zs,
     $                       CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                       CCTK_VARIABLE_REAL,
     $                       in_array_indices(1),
     $                       momentum_int2,
     $                       CCTK_VARIABLE_REAL)


         CALL momentum_integrand3D(origin,Dx,Dy,Dz,x,y,z,
     &        gxx,gxy,gxz,gyy,gyz,gzz,hxx,hxy,hxz,hyy,hyz,hzz,
     &        Extract_temp3d,Psi,Psi_power,conformal_state)
         CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")
         call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                       npoints, 1, 1,
     $                       xs, ys, zs,
     $                       CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                       CCTK_VARIABLE_REAL,
     $                       in_array_indices(1),
     $                       momentum_int3,
     $                       CCTK_VARIABLE_REAL)

      END IF

c     Calculate integrands for spin
c     -----------------------------
      IF (do_spin==1) THEN

         CALL spin_integrand3D(origin,x,y,z,
     &        gxx,gxy,gxz,gyy,gyz,gzz,hxx,hxy,hxz,hyy,hyz,hzz,
     &        Extract_temp3d,Psi,Psi_power,conformal_state)
         CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")
         call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                       npoints, 1, 1,
     $                       xs, ys, zs,
     $                       CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                       CCTK_VARIABLE_REAL,
     $                       in_array_indices(1),
     $                       spin_int1,
     $                       CCTK_VARIABLE_REAL)

         CALL spin_integrand3D(origin,x,y,z,
     &        gxx,gxy,gxz,gyy,gyz,gzz,hxx,hxy,hxz,hyy,hyz,hzz,
     &        Extract_temp3d,Psi,Psi_power,conformal_state)
         CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")
         call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                       npoints, 1, 1,
     $                       xs, ys, zs,
     $                       CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                       CCTK_VARIABLE_REAL,
     $                       in_array_indices(1),
     $                       spin_int2,
     $                       CCTK_VARIABLE_REAL)


         CALL spin_integrand3D(origin,x,y,z,
     &        gxx,gxy,gxz,gyy,gyz,gzz,hxx,hxy,hxz,hyy,hyz,hzz,
     &        Extract_temp3d,Psi,Psi_power,conformal_state)
         CALL CCTK_SyncGroup(ierror,cctkGH,"extract::temps")
         call CCTK_InterpGV (ierror, cctkGH, interp_handle, coord_system_handle,
     $                       npoints, 1, 1,
     $                       xs, ys, zs,
     $                       CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
     $                       CCTK_VARIABLE_REAL,
     $                       in_array_indices(1),
     $                       spin_int3,
     $                       CCTK_VARIABLE_REAL)

      END IF

      END SUBROUTINE D3_to_D2