1 | MODULE volcano_mod |
---|
2 | |
---|
3 | IMPLICIT NONE |
---|
4 | ! saved,protected variables: (local to the core) |
---|
5 | |
---|
6 | integer,save,protected :: nvolc ! number of active volcanoes |
---|
7 | !$OMP THREADPRIVATE(nvolc) |
---|
8 | |
---|
9 | real,save,allocatable,protected :: latlonvals(:,:) ! dimension (nvolc,2) |
---|
10 | ! lat, lon position of each volcano (deg) |
---|
11 | !$OMP THREADPRIVATE(latlonvals) |
---|
12 | |
---|
13 | integer,save,allocatable,protected :: heightvals(:) ! dimension (nvolc) |
---|
14 | ! altitude grid point index |
---|
15 | ! of eruption for each volcano |
---|
16 | !$OMP THREADPRIVATE(heightvals) |
---|
17 | |
---|
18 | real,save,allocatable,protected :: eruptfreqs(:,:) ! dimension (nvolc,3) |
---|
19 | ! for each volcano: interval |
---|
20 | ! between each eruption in days, |
---|
21 | ! day offset from 0, |
---|
22 | ! and duration of eruption in days. |
---|
23 | !$OMP THREADPRIVATE(eruptfreqs) |
---|
24 | |
---|
25 | real,save,allocatable,protected :: volfluxes(:,:) ! dimension (nvolc,nq) |
---|
26 | ! outgassing flux of each tracer for |
---|
27 | ! each volcano, in kg/s |
---|
28 | !$OMP THREADPRIVATE(volfluxes) |
---|
29 | |
---|
30 | integer,save,allocatable,protected :: ivolc(:) ! dimension (nvolc) |
---|
31 | ! horizontal grid index of volcano |
---|
32 | !$OMP THREADPRIVATE(ivolc) |
---|
33 | |
---|
34 | real,save,allocatable,protected :: dtmp(:) ! dimension (nvolc) |
---|
35 | ! time counter for eruption |
---|
36 | !$OMP THREADPRIVATE(dtmp) |
---|
37 | |
---|
38 | logical,save,allocatable,protected :: eruptbool(:) ! dimension (nvolc) |
---|
39 | ! flag to signal whether to erupt volcano or not. |
---|
40 | !$OMP THREADPRIVATE(eruptbool) |
---|
41 | |
---|
42 | CONTAINS |
---|
43 | |
---|
44 | SUBROUTINE volcano(ngrid,nlay,pplev,wu,wv,pt,ssource,nq, & |
---|
45 | massarea,zday,ptimestep) |
---|
46 | |
---|
47 | use time_phylmdz_mod, only: daysec ! day length (s) |
---|
48 | use vertical_layers_mod, only: pseudoalt !atm. layer pseudo-altitude (km) |
---|
49 | |
---|
50 | IMPLICIT NONE |
---|
51 | |
---|
52 | !======================================================================= |
---|
53 | ! Volcanic Activity (updated to run in parallel) |
---|
54 | ! |
---|
55 | ! Description : |
---|
56 | ! Author : Lionel Wilson |
---|
57 | ! Adapted for the LMD/GCM by Laura Kerber, J.-B. Madeleine, Saira Hamid |
---|
58 | ! Adapted to be more generalisable by Ashwin Braude (02.11.2023) |
---|
59 | ! Further optimised for parallelisation by Ehouarn Millour (12.04.2024) |
---|
60 | ! |
---|
61 | ! Reference : |
---|
62 | ! @ARTICLE{2007JVGR..163...83W, |
---|
63 | ! author = {{Wilson}, L. and {Head}, J.~W.}, |
---|
64 | ! title = "{Explosive volcanic eruptions on Mars: Tephra and |
---|
65 | ! accretionary lapilli formation, dispersal and recognition in the |
---|
66 | ! geologic record}", |
---|
67 | ! journal = {Journal of Volcanology and Geothermal Research}, |
---|
68 | ! year = 2007, |
---|
69 | ! month = jun, |
---|
70 | ! volume = 163} |
---|
71 | !======================================================================= |
---|
72 | |
---|
73 | ! Variable statement |
---|
74 | ! __________________________________________________________________ |
---|
75 | |
---|
76 | |
---|
77 | ! Parameters |
---|
78 | ! ---------- |
---|
79 | |
---|
80 | ! PLEASE DEFINE THE LOCATION OF THE ERUPTION BELOW : |
---|
81 | ! ============================================ |
---|
82 | ! Source flux (kg/s) |
---|
83 | ! REAL, parameter :: mmsource = 1E8 !Mastin et al. 2009 |
---|
84 | ! REAL, parameter :: wsource = 1.E6 !1wt% magma water content(Greeley 1987) |
---|
85 | ! REAL, parameter :: h2so4source = 1.E8 |
---|
86 | ! ============================================= |
---|
87 | ! Local variables |
---|
88 | ! --------------- |
---|
89 | |
---|
90 | INTEGER :: l |
---|
91 | INTEGER :: iq ! Tracer identifier |
---|
92 | INTEGER :: ig |
---|
93 | INTEGER :: volci |
---|
94 | |
---|
95 | LOGICAL,SAVE :: firstcall=.true. |
---|
96 | !$OMP THREADPRIVATE(firstcall) |
---|
97 | |
---|
98 | ! Inputs |
---|
99 | ! ------ |
---|
100 | |
---|
101 | INTEGER, intent(in) :: nq ! Number of tracers |
---|
102 | INTEGER, intent(in) :: ngrid ! Number of grid points |
---|
103 | INTEGER,intent(in) :: nlay ! Number of vertical levels |
---|
104 | REAL,intent(in) :: pplev(ngrid,nlay+1) ! Pressure between the layers (Pa) |
---|
105 | REAL,intent(in) :: wv(ngrid,nlay+1) ! wind |
---|
106 | REAL,intent(in) :: wu(ngrid,nlay+1) ! wind |
---|
107 | REAL,intent(in) :: pt(ngrid,nlay+1) ! temp |
---|
108 | REAL,intent(in) :: massarea(ngrid,nlay) |
---|
109 | REAL,intent(in) :: zday ! "Universal time": in sols (and fraction thereof). |
---|
110 | ! since the North. Spring equinox (Ls=0) |
---|
111 | REAL,intent(in) :: ptimestep ! Physics timestep (s). |
---|
112 | |
---|
113 | ! Outputs |
---|
114 | ! ------- |
---|
115 | REAL, intent(out) :: ssource(ngrid,nlay,nq) ! Source tendency (kg.kg-1.s-1) |
---|
116 | |
---|
117 | ! Initialization |
---|
118 | ! __________________________________________________________________ |
---|
119 | |
---|
120 | if (firstcall) then |
---|
121 | ! initialize everything the very first time this routine is called |
---|
122 | call read_volcano(nq,ngrid,nlay) |
---|
123 | !Make sure eruptfreqs values are multiples of timestep |
---|
124 | ! within roundoffs of modulo hence comparison to 1.e-6*(ptimestep/daysec) and not 0 |
---|
125 | do volci=1,nvolc |
---|
126 | ! if((modulo(eruptfreqs(volci,1),(ptimestep/daysec)).ne.0).or.(modulo(eruptfreqs(volci,2),(ptimestep/daysec)).ne.0))then |
---|
127 | if ((modulo(eruptfreqs(volci,1),(ptimestep/daysec)).gt.1.e-6*(ptimestep/daysec)).or. & |
---|
128 | (modulo(eruptfreqs(volci,2),(ptimestep/daysec)).gt.1.e-6*(ptimestep/daysec))) then |
---|
129 | write(*,*) volci,eruptfreqs(volci,1)/(ptimestep/daysec),eruptfreqs(volci,2)/(ptimestep/daysec) |
---|
130 | call abort_physic('volcano','frequency and offset of eruption must both be multiples of timestep',1) |
---|
131 | endif |
---|
132 | if(eruptfreqs(volci,3)-eruptfreqs(volci,1).ge.0.0)then |
---|
133 | write(*,*) volci, eruptfreqs(volci,1), eruptfreqs(volci,3) |
---|
134 | call abort_physic('volcano','duration of eruption must be lower than frequency',1) |
---|
135 | endif |
---|
136 | write(*,*)'volcano: First eruption of volcano ',volci, ' scheduled in ', & |
---|
137 | modulo(zday-eruptfreqs(volci,2),eruptfreqs(volci,1)), 'days' |
---|
138 | enddo |
---|
139 | firstcall=.false. |
---|
140 | endif |
---|
141 | |
---|
142 | ssource(1:ngrid,1:nlay,1:nq)=0 !all arrays=zero since it's a local variable within routine that need to be initialized |
---|
143 | |
---|
144 | do volci=1,nvolc |
---|
145 | |
---|
146 | !If time to start erupting volcano |
---|
147 | if(modulo(nint((zday-eruptfreqs(volci,2))*daysec/ptimestep),nint(eruptfreqs(volci,1)*daysec/ptimestep)).eq.0)then |
---|
148 | |
---|
149 | dtmp(volci) = eruptfreqs(volci,3) |
---|
150 | eruptbool(volci) = .true. |
---|
151 | |
---|
152 | !emit source flux over the duration of the eruption |
---|
153 | if(dtmp(volci).lt.ptimestep/daysec)then |
---|
154 | WRITE(*,*) 'Volcano ',volci,' Time: ', zday, eruptfreqs(volci,:) |
---|
155 | l = heightvals(volci) |
---|
156 | ig=ivolc(volci) |
---|
157 | WRITE(*,*) 'Erupt volcano at ivolc=',ivolc(volci),' height=',pseudoalt(l), ' for ', dtmp(volci), ' days' |
---|
158 | do iq=1,nq |
---|
159 | ssource(ig,l,iq) = ssource(ig,l,iq)+& |
---|
160 | (dtmp(volci)*daysec/ptimestep)*volfluxes(volci,iq)/massarea(ig,l) |
---|
161 | enddo |
---|
162 | eruptbool(volci) = .false. |
---|
163 | endif |
---|
164 | endif |
---|
165 | |
---|
166 | !If eruption lasts for more than a single timestep, need to continue erupting. |
---|
167 | if(eruptbool(volci))then |
---|
168 | WRITE(*,*) 'Volcano ',volci,' Time: ', zday, eruptfreqs(volci,:) |
---|
169 | l = heightvals(volci) |
---|
170 | ig=ivolc(volci) |
---|
171 | WRITE(*,*) 'Erupt volcano at ivolc=',ivolc(volci),' height=',pseudoalt(l), ' for ', dtmp(volci), ' more days' |
---|
172 | |
---|
173 | do iq=1,nq |
---|
174 | ssource(ig,l,iq) = ssource(ig,l,iq)+& |
---|
175 | (MIN(dtmp(volci),ptimestep/daysec)*daysec/ptimestep)*volfluxes(volci,iq)/massarea(ig,l) |
---|
176 | enddo |
---|
177 | |
---|
178 | dtmp(volci) = dtmp(volci) - ptimestep/daysec |
---|
179 | if(dtmp(volci).lt.0.0) eruptbool(volci) = .false. |
---|
180 | endif |
---|
181 | |
---|
182 | enddo ! of do volci=1,nvolc |
---|
183 | |
---|
184 | END SUBROUTINE volcano |
---|
185 | |
---|
186 | !======================================================================= |
---|
187 | |
---|
188 | SUBROUTINE read_volcano(nq,ngrid,nlay) |
---|
189 | ! this routine reads in volcano.def and initializes module variables |
---|
190 | USE mod_phys_lmdz_para, only: is_master, bcast |
---|
191 | use mod_phys_lmdz_transfert_para, only: reduce_sum |
---|
192 | use comcstfi_mod, only: pi |
---|
193 | use tracer_h, only: is_known_tracer, tracer_index |
---|
194 | use mod_grid_phy_lmdz, only: nvertex ! number of grid cell "corners" |
---|
195 | use geometry_mod, only: boundslon, boundslat ! grid cell boundaries (rad) |
---|
196 | use geometry_mod, only: longitude_deg, latitude_deg ! grid cell coordinates (deg) |
---|
197 | use vertical_layers_mod, only: pseudoalt !atm. layer pseudo-altitude (km) |
---|
198 | IMPLICIT NONE |
---|
199 | |
---|
200 | integer, intent(in) :: nq ! number of tracers |
---|
201 | INTEGER, intent(in) :: ngrid ! number of atmospheric columns |
---|
202 | INTEGER, intent(in) :: nlay ! |
---|
203 | real :: zlev(nlay+1) ! inter-layer altitudes (km) |
---|
204 | real, allocatable :: latlonvals_total(:,:) !lat, lon position of each volcano |
---|
205 | integer, allocatable :: heightvals_total(:) ! alt index of eruption |
---|
206 | real, allocatable :: eruptfreqs_total(:,:) ! date and freq of eruption |
---|
207 | real, allocatable :: volfluxes_total(:,:) ! injected tracer rate |
---|
208 | real, allocatable :: dtmp_total(:) |
---|
209 | logical, allocatable :: eruptbool_total(:) |
---|
210 | integer :: ierr ! error return code |
---|
211 | integer :: iq,l,volci,ivoltrac,nvoltrac,blank |
---|
212 | real :: eruptfreq,eruptoffset,eruptdur |
---|
213 | integer :: nvolctrac ! number of outgased tracers |
---|
214 | integer :: ig,i |
---|
215 | integer :: nvolc_total ! total (planet-wide) number of volcanoes |
---|
216 | integer :: nvolc_check |
---|
217 | real :: lat_volc,lon_volc,volflux,alt_volc |
---|
218 | character(len=500) :: voltrac(nq) |
---|
219 | character(len=500) :: tracline ! to read volcano.def line |
---|
220 | real :: boundslon_deg(nvertex) ! cell longitude boundaries (deg) |
---|
221 | real :: boundslat_deg(nvertex) ! cell latitude boundaries (deg) |
---|
222 | real :: minlat, minlon, maxlat, maxlon |
---|
223 | real tmp(nlay+1) |
---|
224 | logical,allocatable :: belongs_here(:) |
---|
225 | integer,allocatable :: volc_index(:) |
---|
226 | character(len=*),parameter :: rname="read_volcano" |
---|
227 | |
---|
228 | ! 1. Only the master reads and processes the volcano.def file |
---|
229 | if (is_master) then |
---|
230 | ! 1.1 Open file and get total number of volcanoes |
---|
231 | open(407, form = 'formatted', status = 'old', & |
---|
232 | file = 'volcano.def', iostat=ierr) |
---|
233 | if (ierr /=0) then |
---|
234 | print*,trim(rname)//': Problem in opening volcano.def' |
---|
235 | stop |
---|
236 | end if |
---|
237 | |
---|
238 | write(*,*)trim(rname)//' Reading file volcano.def' |
---|
239 | DO |
---|
240 | READ(407,'(A)',iostat=ierr) tracline |
---|
241 | IF (ierr==0) THEN |
---|
242 | IF (index(tracline,'#').ne.1) THEN ! Allows arbitary number of comments lines in the header |
---|
243 | !Read in the number of point sources/volcanoes in volcano.def |
---|
244 | READ(tracline,*) nvolc_total |
---|
245 | EXIT |
---|
246 | ENDIF |
---|
247 | ELSE ! If pb, or if reached EOF without having found nqtot |
---|
248 | ! call abort_physic('initracer','Unable to read numbers '// & |
---|
249 | ! 'of tracers in traceur.def',1) |
---|
250 | print*,trim(rname)//": unable to read numbers of tracer in volcano.def" |
---|
251 | stop |
---|
252 | ENDIF |
---|
253 | ENDDO |
---|
254 | endif ! of if (is_master) |
---|
255 | ! tell the world about nvolc_total |
---|
256 | CALL bcast(nvolc_total) |
---|
257 | |
---|
258 | ! 1.2. allocate arrays (all cores) to store all the information |
---|
259 | ALLOCATE(latlonvals_total(nvolc_total,2),stat=ierr) |
---|
260 | if (ierr/=0) then |
---|
261 | write(*,*) trim(rname)//" failed allocation for latlonvals_total()" |
---|
262 | call abort_physic(rname,"failed allocation",1) |
---|
263 | endif |
---|
264 | ALLOCATE(heightvals_total(nvolc_total),stat=ierr) |
---|
265 | if (ierr/=0) then |
---|
266 | write(*,*) trim(rname)//" failed allocation for heightvals_total()" |
---|
267 | call abort_physic(rname,"failed allocation",1) |
---|
268 | endif |
---|
269 | ALLOCATE(eruptfreqs_total(nvolc_total,3),stat=ierr) |
---|
270 | if (ierr/=0) then |
---|
271 | write(*,*) trim(rname)//" failed allocation for eruptfreqs_total()" |
---|
272 | call abort_physic(rname,"failed allocation",1) |
---|
273 | endif |
---|
274 | ALLOCATE(volfluxes_total(nvolc_total,nq),stat=ierr) |
---|
275 | if (ierr/=0) then |
---|
276 | write(*,*) trim(rname)//" failed allocation for volfluxes_total()" |
---|
277 | call abort_physic(rname,"failed allocation",1) |
---|
278 | endif |
---|
279 | ALLOCATE(dtmp_total(nvolc_total),stat=ierr) |
---|
280 | if (ierr/=0) then |
---|
281 | write(*,*) trim(rname)//" failed allocation for dtmp_total()" |
---|
282 | call abort_physic(rname,"failed allocation",1) |
---|
283 | endif |
---|
284 | ALLOCATE(eruptbool_total(nvolc_total),stat=ierr) |
---|
285 | if (ierr/=0) then |
---|
286 | write(*,*) trim(rname)//" failed allocation for eruptbool_total()" |
---|
287 | call abort_physic(rname,"failed allocation",1) |
---|
288 | endif |
---|
289 | ! initialize volfluxes_total |
---|
290 | volfluxes_total(:,:)=0.0 |
---|
291 | |
---|
292 | ! 1.3. continue (master only) to read and process volcano.def |
---|
293 | if (is_master) then |
---|
294 | volci = 0 |
---|
295 | !Iterate for each volcano |
---|
296 | do while (volci.lt.nvolc_total) |
---|
297 | DO |
---|
298 | READ(407,'(A)',iostat=ierr) tracline |
---|
299 | IF (ierr==0) THEN |
---|
300 | IF (index(tracline,'#').ne.1) THEN |
---|
301 | READ(tracline,*) volci |
---|
302 | EXIT |
---|
303 | ENDIF |
---|
304 | ELSE ! If pb, or if reached EOF without having found nqtot |
---|
305 | ! call abort_physic('initracer','Unable to read numbers '// & |
---|
306 | ! 'of tracers in traceur.def',1) |
---|
307 | print*,trim(rname)//": unable to read line in volcano.def" |
---|
308 | stop |
---|
309 | ENDIF |
---|
310 | ENDDO |
---|
311 | write(*,*)trim(rname)//': volcano number:',volci |
---|
312 | READ(407,*,iostat=ierr) lat_volc,lon_volc,alt_volc |
---|
313 | !Convert altitude in pressure coords to altitude in grid coords. |
---|
314 | ! First compute inter-layer pseudo altitudes (m) |
---|
315 | zlev(1)=0. |
---|
316 | do i=2,nlay |
---|
317 | ! inter-layer i is half way between mid-layers i and i-1 |
---|
318 | ! and convert km to m |
---|
319 | zlev(i)=1.e3*(pseudoalt(i)+pseudoalt(i-1))/2. |
---|
320 | enddo |
---|
321 | !extrapolate topmost layer top boundary |
---|
322 | zlev(nlay+1)= 2*zlev(nlay)-zlev(nlay-1) |
---|
323 | ! Then identify eruption layer index |
---|
324 | do i=1,nlay+1 |
---|
325 | tmp(i) = abs(zlev(i)-alt_volc*1000.0) |
---|
326 | enddo |
---|
327 | l = minloc(tmp,1) |
---|
328 | if(zlev(l).gt.alt_volc)then |
---|
329 | l = l-1 |
---|
330 | endif |
---|
331 | write(*,*)trim(rname)//' latitude:',lat_volc |
---|
332 | write(*,*)trim(rname)//' longitude:',lon_volc |
---|
333 | write(*,*)trim(rname)//' altitude snapped to between',zlev(l), ' and ', zlev(l+1), 'm' |
---|
334 | latlonvals_total(volci,1) = lat_volc |
---|
335 | latlonvals_total(volci,2) = lon_volc |
---|
336 | heightvals_total(volci) = l |
---|
337 | !Read in number of active tracers being outgassed |
---|
338 | READ(407,*,iostat=ierr) nvoltrac |
---|
339 | write(*,*)trim(rname)//' number of outgased tracers:',nvoltrac |
---|
340 | !Read in eruption frequency (in multiples of iphysiq) and offset of eruption |
---|
341 | READ(407,'(A)',iostat=ierr) tracline |
---|
342 | READ(tracline,*) eruptfreq, eruptoffset, eruptdur |
---|
343 | write(*,*)trim(rname)//' eruption frequency (sols):',eruptfreq |
---|
344 | write(*,*)trim(rname)//' eruption offset (sols):',eruptoffset |
---|
345 | write(*,*)trim(rname)//' eruption duration (sols):',eruptdur |
---|
346 | eruptfreqs_total(volci,1) = eruptfreq |
---|
347 | eruptfreqs_total(volci,2) = eruptoffset |
---|
348 | eruptfreqs_total(volci,3) = eruptdur |
---|
349 | dtmp_total(volci) = eruptdur |
---|
350 | eruptbool_total(volci) = .false. |
---|
351 | |
---|
352 | do ivoltrac=1,nvoltrac |
---|
353 | |
---|
354 | !Read in tracer name and check it's in traceur.def |
---|
355 | read(407,'(A)') tracline |
---|
356 | blank = index(tracline,' ') ! Find position of 1st blank in tracline |
---|
357 | voltrac(ivoltrac) = tracline(1:blank) |
---|
358 | ! check it is a valid tracer: |
---|
359 | if (is_known_tracer(voltrac(ivoltrac))) then |
---|
360 | !Read in flux of tracer out of volcano (kg/s) |
---|
361 | READ(407,'(A)',iostat=ierr) tracline |
---|
362 | READ(tracline,*) volflux |
---|
363 | iq=tracer_index(voltrac(ivoltrac)) |
---|
364 | volfluxes_total(volci,iq) = volflux |
---|
365 | else |
---|
366 | write(*,*) trim(rname)//": problem with ",trim(voltrac(ivoltrac)) |
---|
367 | call abort_physic(rname, 'tracer '//& |
---|
368 | trim(voltrac(ivoltrac))//' not found', 1) |
---|
369 | endif |
---|
370 | |
---|
371 | enddo ! of do ivoltrac=1,nvoltrac |
---|
372 | |
---|
373 | enddo ! of do while (volci.lt.nvolc_total) |
---|
374 | |
---|
375 | close(407) |
---|
376 | endif ! of if (is_master) |
---|
377 | |
---|
378 | ! 2. Share data with all cores |
---|
379 | call bcast(latlonvals_total) |
---|
380 | call bcast(heightvals_total) |
---|
381 | call bcast(eruptfreqs_total) |
---|
382 | call bcast(volfluxes_total) |
---|
383 | call bcast(dtmp_total) |
---|
384 | call bcast(eruptbool_total) |
---|
385 | |
---|
386 | ! 3. Each core extracts and stores data relevant to its domain only |
---|
387 | ! 3.1. find out how many volcanoes for this core |
---|
388 | nvolc=0 |
---|
389 | allocate(belongs_here(nvolc_total),stat=ierr) |
---|
390 | if (ierr/=0) then |
---|
391 | write(*,*) trim(rname)//" failed allocation for belong_here()" |
---|
392 | call abort_physic(rname,"failled allocation",1) |
---|
393 | endif |
---|
394 | allocate(volc_index(nvolc_total),stat=ierr) |
---|
395 | if (ierr/=0) then |
---|
396 | write(*,*) trim(rname)//" failed allocation for volc_index()" |
---|
397 | call abort_physic(rname,"failled allocation",1) |
---|
398 | endif |
---|
399 | |
---|
400 | do volci=1,nvolc_total |
---|
401 | belongs_here(volci)=.false. |
---|
402 | volc_index(volci)=0 |
---|
403 | |
---|
404 | lat_volc=latlonvals_total(volci,1) |
---|
405 | lon_volc=latlonvals_total(volci,2) |
---|
406 | |
---|
407 | ! find out if this volcano is in the ig cell |
---|
408 | do ig=1,ngrid |
---|
409 | boundslon_deg(:)=boundslon(ig,:)/pi*180. |
---|
410 | boundslat_deg(:)=boundslat(ig,:)/pi*180. |
---|
411 | |
---|
412 | minlon=minval(boundslon_deg) |
---|
413 | maxlon=maxval(boundslon_deg) |
---|
414 | minlat=minval(boundslat_deg) |
---|
415 | maxlat=maxval(boundslat_deg) |
---|
416 | |
---|
417 | if ((minlat<=lat_volc).and.(lat_volc<=maxlat).and. & |
---|
418 | (((minlon<=lon_volc).and.(lon_volc<=maxlon)).or. & |
---|
419 | ((minlon+360.<=lon_volc).and.(lon_volc<=maxlon+360.)).or. & |
---|
420 | ((minlon-360.<=lon_volc).and.(lon_volc<=maxlon-360.)))) then |
---|
421 | nvolc=nvolc+1 |
---|
422 | belongs_here(volci)=.true. |
---|
423 | volc_index(volci)=ig |
---|
424 | ! the job is done move on to next volcano |
---|
425 | exit |
---|
426 | endif |
---|
427 | enddo ! of do ig=1,ngrid |
---|
428 | enddo ! of volci=1,nvolc_total |
---|
429 | |
---|
430 | ! 3.2. allocate arrays, now that nvolc is known |
---|
431 | if (nvolc>0) then |
---|
432 | allocate(latlonvals(nvolc,2),stat=ierr) |
---|
433 | if (ierr/=0) then |
---|
434 | write(*,*) trim(rname)//" failed allocation for latlonvals()" |
---|
435 | call abort_physic(rname,"failed allocation",1) |
---|
436 | endif |
---|
437 | allocate(heightvals(nvolc),stat=ierr) |
---|
438 | if (ierr/=0) then |
---|
439 | write(*,*) trim(rname)//" failed allocation for heightvals()" |
---|
440 | call abort_physic(rname,"failed allocation",1) |
---|
441 | endif |
---|
442 | allocate(eruptfreqs(nvolc,3),stat=ierr) |
---|
443 | if (ierr/=0) then |
---|
444 | write(*,*) trim(rname)//" failed allocation for eruptfreqs()" |
---|
445 | call abort_physic(rname,"failed allocation",1) |
---|
446 | endif |
---|
447 | allocate(volfluxes(nvolc,nq),stat=ierr) |
---|
448 | if (ierr/=0) then |
---|
449 | write(*,*) trim(rname)//" failed allocation for volfluxes()" |
---|
450 | call abort_physic(rname,"failed allocation",1) |
---|
451 | endif |
---|
452 | allocate(ivolc(nvolc),stat=ierr) |
---|
453 | if (ierr/=0) then |
---|
454 | write(*,*) trim(rname)//" failed allocation for ivolc()" |
---|
455 | call abort_physic(rname,"failed allocation",1) |
---|
456 | endif |
---|
457 | allocate(dtmp(nvolc),stat=ierr) |
---|
458 | if (ierr/=0) then |
---|
459 | write(*,*) trim(rname)//" failed allocation for dtmp()" |
---|
460 | call abort_physic(rname,"failed allocation",1) |
---|
461 | endif |
---|
462 | allocate(eruptbool(nvolc),stat=ierr) |
---|
463 | if (ierr/=0) then |
---|
464 | write(*,*) trim(rname)//" failed allocation for eruptbool()" |
---|
465 | call abort_physic(rname,"failed allocation",1) |
---|
466 | endif |
---|
467 | endif ! of if (nvolc>0) |
---|
468 | |
---|
469 | ! 3.3. copy over global data to local arrays |
---|
470 | do volci=1,nvolc |
---|
471 | do i=1,nvolc_total |
---|
472 | if (belongs_here(i)) then |
---|
473 | ! copy data to local arrays |
---|
474 | latlonvals(volci,1:2)=latlonvals_total(i,1:2) |
---|
475 | heightvals(volci)=heightvals_total(i) |
---|
476 | eruptfreqs(volci,1:3)=eruptfreqs_total(i,1:3) |
---|
477 | volfluxes(volci,1:nq)=volfluxes_total(i,1:nq) |
---|
478 | ivolc(volci)=volc_index(i) |
---|
479 | dtmp(volci) = dtmp_total(i) |
---|
480 | eruptbool(volci) = eruptbool_total(i) |
---|
481 | |
---|
482 | ! the job is done, move on to next one |
---|
483 | belongs_here(i)=.false. |
---|
484 | exit |
---|
485 | endif ! of if (belongs_here(i)) |
---|
486 | enddo ! of do i=1,nvolc_total |
---|
487 | enddo ! of do volci=1,nvolc |
---|
488 | |
---|
489 | ! 4. Sanity check |
---|
490 | ! verify that the sum of volcanoes over all cores is indeed nvolc_total |
---|
491 | if (is_master) nvolc_check=0 |
---|
492 | call bcast(nvolc_check) |
---|
493 | call reduce_sum(nvolc,nvolc_check) |
---|
494 | call bcast(nvolc_check) |
---|
495 | if (nvolc_check/=nvolc_total) then |
---|
496 | write(*,*) trim(rname)//" error. Sum of volcanoes over all cores ", & |
---|
497 | nvolc_check," differs from nvolc_total=",nvolc_total |
---|
498 | call abort_physic(rname," error. Check and fix volcano.F90",1) |
---|
499 | endif |
---|
500 | |
---|
501 | ! give a summary, for each core: |
---|
502 | write(*,*) trim(rname)//": nvolc=",nvolc," /",nvolc_total |
---|
503 | do volci=1,nvolc |
---|
504 | write(*,*) trim(rname)//": latlonvals(1:2)=",latlonvals(volci,1:2) |
---|
505 | write(*,*) trim(rname)//": heightvals=",heightvals(volci) |
---|
506 | write(*,*) trim(rname)//": eruptfreqs(1:3)=",eruptfreqs(volci,1:3) |
---|
507 | write(*,*) trim(rname)//": volfluxes(1:nq)=",volfluxes(volci,1:nq) |
---|
508 | write(*,*) trim(rname)//": ivolc=",ivolc(volci) |
---|
509 | write(*,*) trim(rname)//": longitude(ivolc)=",longitude_deg(ivolc(volci)) |
---|
510 | write(*,*) trim(rname)//": latitude(ivolc)=",latitude_deg(ivolc(volci)) |
---|
511 | enddo |
---|
512 | |
---|
513 | END SUBROUTINE read_volcano |
---|
514 | |
---|
515 | END MODULE volcano_mod |
---|