print_starwall_response.f90 12.1 KB
Newer Older
Serhiy Mochalskyy's avatar
Serhiy Mochalskyy committed
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
subroutine print_starwall_response

use icontr
use contr_su
use solv 
use tri_w
use tri_p
use coil2d

use sca
use mpi_v
use time
use resistive
!-----------------------------------------------------------------------
implicit none
include "mpif.h"


integer  :: i,j,k,i_loc,j_loc
real     :: num=0.
real,dimension(:,:),allocatable ::  a_ye_print(:,:),a_ey_print(:,:),a_ee_print(:,:)
real,dimension(:,:),allocatable ::  s_ww_print(:,:),s_ww_inv_print(:,:)
!-----------------------------------------------------------------------

! (Determine positions of wall triangle nodes)
if(rank==0) then
   allocate( xyzpot_w(npot_w,3) )
   do i = 1, ntri_w
      do j = 1, 3
         xyzpot_w(jpot_w(i,j),:) = (/ xw(i,j), yw(i,j), zw(i,j) /)
      end do
   end do
endif

if (format_type == 'formatted' ) then
  if(rank==0) then 
       open(60, file='starwall-response.dat', form=format_type, status='replace', action='write')
       130 format(a)
       131 format('#@intparam  ',a24,99i12)
       132 format('#@array     ',a24,i12,a24,99i12)
       133 format(4es24.16)
       134 format(8i12)

       write(60,130) '#@content STARWALL VACUUM RESPONSE DATA FOR THE JOREK CODE'
       write(60,131) 'file_version', 1
       write(60,131) 'n_bnd', N_bnd
       write(60,131) 'nd_bez', nd_bez
       write(60,131) 'ncoil', ncoil
       write(60,131) 'npot_w', npot_w
       write(60,131) 'n_w', n_w
       write(60,131) 'ntri_w', ntri_w
       write(60,131) 'n_tor', n_tor_jorek
       write(60,132) 'i_tor', 1, 'int', n_tor_jorek, 0
       write(60,134) i_tor_jorek(1:n_tor_jorek)
       write(60,132) 'yy', 1, 'float', n_w, 0
       write(60,133) 1./gamma(:)
       write(60,132) 'ye', 2, 'float', n_w, nd_bez

       allocate(a_ye_print(n_w,nd_bez),stat=ier)
           IF (IER /= 0) THEN
               WRITE (*,*) "output.f90, can not allocate a_ye_print  MYPROC_NUM=",MYPNUM
               STOP
           END IF

   endif

!========================================================================================= 
   CALL DESCINIT(DESCA, n_w,nd_bez, NB, NB, 0, 0, CONTEXT, LDA_ye, INFO_A )
   if(INFO_A .NE. 0) then
          write(6,*) "Something is wrong in output  CALL DESCINIT DESCA, INFO_A=",INFO_A
          stop
   endif


   DO i_loc = 1,n_w
       DO j_loc = 1,nd_bez
            CALL pdelget('A','D',num, a_ye_loc,i_loc,j_loc,DESCA)
            if(rank==0) a_ye_print(i_loc,j_loc)=num
       END DO
   END DO
   

  if(rank==0) then
     write(60,133) a_ye_print(:,:)
  
     deallocate (a_ye_print)
     allocate(a_ey_print(nd_bez,n_w),stat=ier)
           IF (IER /= 0) THEN
               WRITE (*,*) "output.f90, can not allocate a_ey_print  MYPROC_NUM=",MYPNUM
               STOP
           END IF
  endif
!========================================================================================= 
   CALL DESCINIT(DESCA, nd_bez,n_w, NB, NB, 0, 0, CONTEXT, LDA_ey, INFO_A )
   if(INFO_A .NE. 0) then
          write(6,*) "Something is wrong in output  CALL DESCINIT DESCA, INFO_A=",INFO_A
          stop
   endif

   DO i_loc = 1,nd_bez
       DO j_loc = 1,n_w
            CALL pdelget('A','D',num, a_ey_loc,i_loc,j_loc,DESCA)
            if(rank==0) a_ey_print(i_loc,j_loc)=num
       END DO
   END DO

  if(rank==0) then
     write(60,132) 'ey', 2, 'float', nd_bez, n_w
     write(60,133) a_ey_print(:,:)
     deallocate (a_ey_print)   
     allocate(a_ee_print(nd_bez,nd_bez),stat=ier)
           IF (IER /= 0) THEN
               WRITE (*,*) "output.f90, can not allocate a_ee_print  MYPROC_NUM=",MYPNUM
               STOP
           END IF

  endif
!========================================================================================= 

   CALL DESCINIT(DESCA, nd_bez,nd_bez, NB, NB, 0, 0, CONTEXT, LDA_ee, INFO_A )
   if(INFO_A .NE. 0) then
          write(6,*) "Something is wrong in output  CALL DESCINIT DESCA, INFO_A=",INFO_A
          stop
   endif

   DO i_loc = 1,nd_bez
       DO j_loc = 1,nd_bez
            CALL pdelget('A','D',num, a_ee_loc,i_loc,j_loc,DESCA)
            if(rank==0) a_ee_print(i_loc,j_loc)=num
       END DO
   END DO


  if(rank==0) then
     write(60,132) 'ee', 2, 'float', nd_bez, nd_bez
     write(60,133) a_ee_print(:,:)
     deallocate (a_ee_print)
     allocate(s_ww_print(n_w,n_w),stat=ier)
           IF (IER /= 0) THEN
               WRITE (*,*) "output.f90, can not allocate s_ww_print  MYPROC_NUM=",MYPNUM
               STOP
           END IF

  endif
!========================================================================================= 

   CALL DESCINIT(DESCA, n_w,n_w, NB, NB, 0, 0, CONTEXT, LDA_sww, INFO_A )
   if(INFO_A .NE. 0) then
          write(6,*) "Something is wrong in output  CALL DESCINIT DESCA, INFO_A=",INFO_A
          stop
   endif

   DO i_loc = 1,n_w
       DO j_loc = 1,n_w
            CALL pdelget('A','D',num, s_ww_loc,i_loc,j_loc,DESCA)
            if(rank==0) s_ww_print(i_loc,j_loc)=num
       END DO
   END DO

  if(rank==0) then
     write(60,132) 's_ww', 2, 'float', n_w, n_w
     write(60,133) s_ww_print(:,:)
     deallocate (s_ww_print)
     allocate(s_ww_inv_print(n_w,n_w),stat=ier)
           IF (IER /= 0) THEN
               WRITE (*,*) "output.f90, can not allocate s_ww_inv_print  MYPROC_NUM=",MYPNUM
               STOP
           END IF

  endif
!========================================================================================= 

 CALL DESCINIT(DESCA, n_w,n_w, NB, NB, 0, 0, CONTEXT, LDA_s_ww_inv, INFO_A )
   if(INFO_A .NE. 0) then
          write(6,*) "Something is wrong in output  CALL DESCINIT DESCA, INFO_A=",INFO_A
          stop
   endif

   DO i_loc = 1,n_w
       DO j_loc = 1,n_w
            CALL pdelget('A','D',num, s_ww_inv_loc,i_loc,j_loc,DESCA)
            if(rank==0) s_ww_inv_print(i_loc,j_loc)=num
       END DO
   END DO

  if(rank==0) then
     write(60,132) 's_ww_inv', 2, 'float', n_w, n_w
     write(60,133) s_ww_inv_print(:,:)
     deallocate (s_ww_inv_print)
  endif
!========================================================================================= 

  if(rank==0) then
     write(60,132) 'xyzpot_w', 2, 'float', npot_w, 3
     write(60,133) xyzpot_w(:,:)
     write(60,132) 'jpot_w', 2, 'int', ntri_w, 3
     write(60,134) jpot_w(:,:)
  endif

else  !if (format_type == 'formatted' ) then


   if(rank==0) then
  
         open(60, file='starwall-response.dat', form=format_type, status='replace', action='write')
         char512='#@content STARWALL VACUUM RESPONSE DATA FOR THE JOREK CODE'
         write(60) 42, 42.d0 !### for an elementary check in JOREK
         write(60) char512
         write(60) '#@intparam  ', 'file_version            ', 1
         write(60) '#@intparam  ', 'n_bnd                   ', N_bnd
         write(60) '#@intparam  ', 'nd_bez                  ', nd_bez
         write(60) '#@intparam  ', 'ncoil                   ', ncoil
         write(60) '#@intparam  ', 'npot_w                  ', npot_w
         write(60) '#@intparam  ', 'n_w                     ', n_w
         write(60) '#@intparam  ', 'ntri_w                  ', ntri_w
         write(60) '#@intparam  ', 'n_tor                   ', n_tor_jorek
         write(60) '#@array     ', 'i_tor                   ', 1, 'int                     ', n_tor_jorek, 0
         write(60) i_tor_jorek(1:n_tor_jorek)
         write(60) '#@array     ', 'yy                      ', 1, 'float                   ', n_w, 0
         write(60) 1.d0/gamma(:)
         write(60) '#@array     ', 'ye                      ', 2, 'float                   ', n_w, nd_bez

         allocate(a_ye_print(n_w,nd_bez),stat=ier)
         IF (IER /= 0) THEN
               WRITE (*,*) "output.f90, can not allocate a_ye_print  MYPROC_NUM=",MYPNUM
               STOP
         END IF
   
    endif
    !========================================================================================= 
    CALL DESCINIT(DESCA, n_w,nd_bez, NB, NB, 0, 0, CONTEXT, LDA_ye, INFO_A )
    if(INFO_A .NE. 0) then
          write(6,*) "Something is wrong in output  CALL DESCINIT DESCA, INFO_A=",INFO_A
          stop
    endif


    DO i_loc = 1,n_w
       DO j_loc = 1,nd_bez
            CALL pdelget('A','D',num, a_ye_loc,i_loc,j_loc,DESCA)
            if(rank==0) a_ye_print(i_loc,j_loc)=num
       END DO
    END DO


     if(rank==0) then
        write(60) a_ye_print(:,:)
        deallocate (a_ye_print)
        allocate(a_ey_print(nd_bez,n_w),stat=ier)
           IF (IER /= 0) THEN
               WRITE (*,*) "output.f90, can not allocate a_ey_print  MYPROC_NUM=",MYPNUM
               STOP
           END IF
    endif
   !========================================================================================= 
   CALL DESCINIT(DESCA, nd_bez,n_w, NB, NB, 0, 0, CONTEXT, LDA_ey, INFO_A )
   if(INFO_A .NE. 0) then
          write(6,*) "Something is wrong in output  CALL DESCINIT DESCA, INFO_A=",INFO_A
          stop
   endif

   DO i_loc = 1,nd_bez
       DO j_loc = 1,n_w
            CALL pdelget('A','D',num, a_ey_loc,i_loc,j_loc,DESCA)
            if(rank==0) a_ey_print(i_loc,j_loc)=num
       END DO
   END DO

  if(rank==0) then
     write(60) '#@array     ', 'ey                      ', 2, 'float                   ', nd_bez, n_w
     write(60) a_ey_print(:,:)
     deallocate (a_ey_print)
     allocate(a_ee_print(nd_bez,nd_bez),stat=ier)
           IF (IER /= 0) THEN
               WRITE (*,*) "output.f90, can not allocate a_ee_print  MYPROC_NUM=",MYPNUM
               STOP
           END IF

  endif
!========================================================================================= 
  

 CALL DESCINIT(DESCA, nd_bez,nd_bez, NB, NB, 0, 0, CONTEXT, LDA_ee, INFO_A )
   if(INFO_A .NE. 0) then
          write(6,*) "Something is wrong in output  CALL DESCINIT DESCA, INFO_A=",INFO_A
          stop
   endif

   DO i_loc = 1,nd_bez
       DO j_loc = 1,nd_bez
            CALL pdelget('A','D',num, a_ee_loc,i_loc,j_loc,DESCA)
            if(rank==0) a_ee_print(i_loc,j_loc)=num
       END DO
   END DO


  if(rank==0) then
     write(60) '#@array     ', 'ee                      ', 2, 'float                   ', nd_bez, nd_bez
     write(60) a_ee_print(:,:)
     deallocate (a_ee_print)
     allocate(s_ww_print(n_w,n_w),stat=ier)
           IF (IER /= 0) THEN
               WRITE (*,*) "output.f90, can not allocate s_ww_print  MYPROC_NUM=",MYPNUM
               STOP
           END IF

  endif
!========================================================================================= 

   CALL DESCINIT(DESCA, n_w,n_w, NB, NB, 0, 0, CONTEXT, LDA_sww, INFO_A )
   if(INFO_A .NE. 0) then
          write(6,*) "Something is wrong in output  CALL DESCINIT DESCA, INFO_A=",INFO_A
          stop
   endif

   DO i_loc = 1,n_w
       DO j_loc = 1,n_w
            CALL pdelget('A','D',num, s_ww_loc,i_loc,j_loc,DESCA)
            if(rank==0) s_ww_print(i_loc,j_loc)=num
       END DO
   END DO

  if(rank==0) then
     write(60) '#@array     ', 's_ww                    ', 2, 'float                   ', n_w, n_w
     write(60) s_ww_print(:,:)
     deallocate (s_ww_print)
     allocate(s_ww_inv_print(n_w,n_w),stat=ier)
           IF (IER /= 0) THEN
               WRITE (*,*) "output.f90, can not allocate s_ww_inv_print  MYPROC_NUM=",MYPNUM
               STOP
           END IF

  endif
!========================================================================================= 
CALL DESCINIT(DESCA, n_w,n_w, NB, NB, 0, 0, CONTEXT, LDA_s_ww_inv, INFO_A )
   if(INFO_A .NE. 0) then
          write(6,*) "Something is wrong in output  CALL DESCINIT DESCA, INFO_A=",INFO_A
          stop
   endif

   DO i_loc = 1,n_w
       DO j_loc = 1,n_w
            CALL pdelget('A','D',num, s_ww_inv_loc,i_loc,j_loc,DESCA)
            if(rank==0) s_ww_inv_print(i_loc,j_loc)=num
       END DO
   END DO

  if(rank==0) then
     write(60) '#@array     ', 's_ww_inv                ', 2, 'float                   ', n_w, n_w
     write(60) s_ww_inv_print(:,:)
     deallocate (s_ww_inv_print)
  endif
!========================================================================================= 

     if(rank==0) then
        write(60) '#@array     ', 'xyzpot_w                ', 2, 'float                   ', npot_w, 3
        write(60) xyzpot_w(:,:)
        write(60) '#@array     ', 'jpot_w                  ', 2, 'int                     ', ntri_w, 3
        write(60) jpot_w(:,:)
     endif

end if !if (format_type == 'formatted' ) then

close(60)

end subroutine print_starwall_response