This repository has been archived by the owner on Jul 29, 2024. It is now read-only.
-
Notifications
You must be signed in to change notification settings - Fork 69
/
VisualizationOutput.pyx
182 lines (142 loc) · 6.16 KB
/
VisualizationOutput.pyx
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
import numpy as np
from mpi4py import MPI
cimport numpy as np
import os
import cython
try:
import cPickle as pickle
except:
import pickle as pickle # for Python 3 users
cimport Grid
cimport ReferenceState
cimport DiagnosticVariables
cimport PrognosticVariables
cimport ParallelMPI
cdef class VisualizationOutput:
def __init__(self, dict namelist, ParallelMPI.ParallelMPI Pa):
self.uuid = str(namelist['meta']['uuid'])
try:
outpath = str(os.path.join(str(namelist['output']['output_root'])
+ 'Output.' + str(namelist['meta']['simname']) + '.' + self.uuid[-5:]))
self.vis_path = os.path.join(outpath, 'Visualization')
except:
self.vis_path = './Visualization.' + self.uuid[-5:]
if Pa.rank == 0:
try:
os.mkdir(outpath)
except:
pass
try:
os.mkdir(self.vis_path)
except:
pass
try:
self.frequency = namelist['visualization']['frequency']
except:
self.frequency = 1e6
return
cpdef initialize(self):
self.last_vis_time = 0.0
return
@cython.boundscheck(False)
@cython.wraparound(False)
@cython.cdivision(True)
cpdef write(self, Grid.Grid Gr, ReferenceState.ReferenceState RS,
PrognosticVariables.PrognosticVariables PV, DiagnosticVariables.DiagnosticVariables DV,
ParallelMPI.ParallelMPI Pa):
cdef:
double [:,:] local_lwp = np.zeros((Gr.dims.n[0], Gr.dims.n[1]), dtype=np.double, order='c')
double [:,:] reduced_lwp = np.zeros((Gr.dims.n[0], Gr.dims.n[1]), dtype=np.double, order='c')
Py_ssize_t i,j,k,ijk
Py_ssize_t imin = Gr.dims.gw
Py_ssize_t jmin = Gr.dims.gw
Py_ssize_t kmin = Gr.dims.gw
Py_ssize_t imax = Gr.dims.nlg[0] - Gr.dims.gw
Py_ssize_t jmax = Gr.dims.nlg[1] - Gr.dims.gw
Py_ssize_t kmax = Gr.dims.nlg[2] - Gr.dims.gw
Py_ssize_t istride = Gr.dims.nlg[1] * Gr.dims.nlg[2]
Py_ssize_t jstride = Gr.dims.nlg[2]
Py_ssize_t ishift, jshift
Py_ssize_t global_shift_i = Gr.dims.indx_lo[0]
Py_ssize_t global_shift_j = Gr.dims.indx_lo[1]
Py_ssize_t global_shift_k = Gr.dims.indx_lo[2]
Py_ssize_t var_shift
Py_ssize_t i2d, j2d
double dz = Gr.dims.dx[2]
dict out_dict = {}
comm = MPI.COMM_WORLD
try:
var_shift = DV.get_varshift(Gr, 'ql')
with nogil:
for i in xrange(imin, imax):
ishift = i * istride
for j in xrange(jmin, jmax):
jshift = j * jstride
for k in xrange(kmin, kmax):
ijk = ishift + jshift + k
i2d = global_shift_i + i - Gr.dims.gw
j2d = global_shift_j + j - Gr.dims.gw
local_lwp[i2d, j2d] += (RS.rho0[k] * DV.values[var_shift + ijk] * dz)
comm.Reduce(local_lwp, reduced_lwp, op=MPI.SUM)
del local_lwp
if Pa.rank == 0:
out_dict['lwp'] = np.array(reduced_lwp,dtype=np.double)
del reduced_lwp
except:
Pa.root_print('Trouble Writing LWP')
#Write output of Prognostic Variables and Diagnostic Variables
cdef:
double [:,:] local_var
double [:,:] reduced_var
list pv_vars = ['qt', 's', 'w']
list dv_vars = ['ql', 'diffusivity']
for var in pv_vars:
local_var = np.zeros((Gr.dims.n[1], Gr.dims.n[2]), dtype=np.double, order='c')
reduced_var = np.zeros((Gr.dims.n[1], Gr.dims.n[2]), dtype=np.double, order='c')
try:
var_shift = PV.get_varshift(Gr, var)
with nogil:
if global_shift_i == 0:
i = 0
ishift = i * istride
for j in xrange(jmin, jmax):
jshift = j * jstride
for k in xrange(kmin, kmax):
ijk = ishift + jshift + k
j2d = global_shift_j + j - Gr.dims.gw
k2d = global_shift_k + k - Gr.dims.gw
local_var[j2d, k2d] = PV.values[var_shift + ijk]
comm.Reduce(local_var, reduced_var, op=MPI.SUM)
del local_var
if Pa.rank == 0:
out_dict[var] = np.array(reduced_var, dtype=np.double)
del reduced_var
except:
Pa.root_print('Trouble Writing ' + var)
for var in dv_vars:
local_var = np.zeros((Gr.dims.n[1], Gr.dims.n[2]), dtype=np.double, order='c')
reduced_var = np.zeros((Gr.dims.n[1], Gr.dims.n[2]), dtype=np.double, order='c')
try:
var_shift = DV.get_varshift(Gr, var)
with nogil:
if global_shift_i == 0:
i = 0
ishift = i * istride
for j in xrange(jmin, jmax):
jshift = j * jstride
for k in xrange(kmin, kmax):
ijk = ishift + jshift + k
j2d = global_shift_j + j - Gr.dims.gw
k2d = global_shift_k + k - Gr.dims.gw
local_var[j2d, k2d] = DV.values[var_shift + ijk]
comm.Reduce(local_var, reduced_var, op=MPI.SUM)
del local_var
if Pa.rank == 0:
out_dict[var] = np.array(reduced_var, dtype=np.double)
del reduced_var
except:
Pa.root_print('Trouble Writing ' + var)
if Pa.rank == 0:
with open(self.vis_path+ '/' + str(10000000 + np.int(self.last_vis_time)) + '.pkl', 'wb') as f:
pickle.dump(out_dict, f, protocol=2)
return