regression.F90 coverage: 85.71 %func 80.16 %block
1) module Regression_module
2)
3) use Output_Aux_module
4)
5) use PFLOTRAN_Constants_module
6)
7) implicit none
8)
9) private
10)
11) #include "petsc/finclude/petscsys.h"
12) #include "petsc/finclude/petscvec.h"
13) #include "petsc/finclude/petscvec.h90"
14)
15) type, public :: regression_type
16) type(regression_variable_type), pointer :: variable_list
17) PetscInt, pointer :: natural_cell_ids(:)
18) PetscInt :: num_cells_per_process
19) PetscInt, pointer :: cells_per_process_natural_ids(:)
20) Vec :: natural_cell_id_vec
21) Vec :: cells_per_process_vec
22) VecScatter :: scatter_natural_cell_id_gtos
23) VecScatter :: scatter_cells_per_process_gtos
24) type(regression_type), pointer :: next
25) end type regression_type
26)
27) type, public :: regression_variable_type
28) character(len=MAXSTRINGLENGTH) :: name
29) type(regression_variable_type), pointer :: next
30) end type regression_variable_type
31)
32) public :: RegressionRead, &
33) RegressionCreateMapping, &
34) RegressionOutput, &
35) RegressionDestroy
36)
37) contains
38)
39) ! ************************************************************************** !
40)
41) function RegressionCreate()
42) !
43) ! Creates a regression object
44) !
45) ! Author: Glenn Hammond
46) ! Date: 10/11/12
47) !
48)
49) implicit none
50)
51) type(regression_type), pointer :: RegressionCreate
52)
53) type(regression_type), pointer :: regression
54)
55) allocate(regression)
56) nullify(regression%variable_list)
57) nullify(regression%natural_cell_ids)
58) regression%num_cells_per_process = 0
59) nullify(regression%cells_per_process_natural_ids)
60) regression%natural_cell_id_vec = 0
61) regression%cells_per_process_vec = 0
62) regression%scatter_natural_cell_id_gtos = 0
63) regression%scatter_cells_per_process_gtos = 0
64) nullify(regression%next)
65) RegressionCreate => regression
66)
67) end function RegressionCreate
68)
69) ! ************************************************************************** !
70)
71) function RegressionVariableCreate()
72) !
73) ! Creates a regression variable object
74) !
75) ! Author: Glenn Hammond
76) ! Date: 10/11/12
77) !
78)
79) implicit none
80)
81) type(regression_variable_type), pointer :: RegressionVariableCreate
82)
83) type(regression_variable_type), pointer :: regression_variable
84)
85) allocate(regression_variable)
86) regression_variable%name = ''
87) nullify(regression_variable%next)
88) RegressionVariableCreate => regression_variable
89)
90) end function RegressionVariableCreate
91)
92) ! ************************************************************************** !
93)
94) subroutine RegressionRead(regression,input,option)
95) !
96) ! Reads in contents of a regression card
97) !
98) ! Author: Glenn Hammond
99) ! Date: 10/11/12
100) !
101)
102) use Option_module
103) use Input_Aux_module
104) use String_module
105) use Utility_module
106)
107) implicit none
108)
109) type(regression_type), pointer :: regression
110) type(input_type), pointer :: input
111) type(option_type) :: option
112)
113) character(len=MAXWORDLENGTH) :: keyword, word
114) type(regression_variable_type), pointer :: cur_variable, new_variable
115) PetscInt :: count, max_cells
116) PetscInt, pointer :: int_array(:)
117) PetscErrorCode :: ierr
118)
119) regression => RegressionCreate()
120)
121) input%ierr = 0
122) do
123)
124) call InputReadPflotranString(input,option)
125)
126) if (InputCheckExit(input,option)) exit
127)
128) call InputReadWord(input,option,keyword,PETSC_TRUE)
129) call InputErrorMsg(input,option,'keyword','REGRESSION')
130) call StringToUpper(keyword)
131)
132) select case(trim(keyword))
133)
134) case('VARIABLES')
135) count = 0
136) do
137) call InputReadPflotranString(input,option)
138) if (InputCheckExit(input,option)) exit
139)
140) call InputReadWord(input,option,word,PETSC_TRUE)
141) call InputErrorMsg(input,option,'variable','REGRESSION,VARIABLES')
142) call StringToUpper(word)
143) new_variable => RegressionVariableCreate()
144) new_variable%name = word
145) if (.not.associated(regression%variable_list)) then
146) regression%variable_list => new_variable
147) else
148) cur_variable%next => new_variable
149) endif
150) cur_variable => new_variable
151) enddo
152) case('CELLS')
153) max_cells = 100
154) allocate(int_array(max_cells))
155) count = 0
156) do
157) call InputReadPflotranString(input,option)
158) if (InputCheckExit(input,option)) exit
159)
160) count = count + 1
161) if (count > max_cells) then
162) call reallocateIntArray(int_array,max_cells)
163) endif
164) call InputReadInt(input,option,int_array(count))
165) call InputErrorMsg(input,option,'natural cell id','REGRESSION,CELLS')
166) enddo
167) allocate(regression%natural_cell_ids(count))
168) regression%natural_cell_ids = int_array(1:count)
169) call PetscSortInt(count,regression%natural_cell_ids, &
170) ierr);CHKERRQ(ierr)
171) deallocate(int_array)
172) case('CELLS_PER_PROCESS')
173) call InputReadInt(input,option,regression%num_cells_per_process)
174) call InputErrorMsg(input,option,'num cells per process','REGRESSION')
175) case default
176) call InputKeywordUnrecognized(keyword,'REGRESSION',option)
177) end select
178)
179) enddo
180)
181) end subroutine RegressionRead
182)
183) ! ************************************************************************** !
184)
185) subroutine RegressionCreateMapping(regression,realization)
186) !
187) ! Creates mapping between a natural mpi vec and a
188) ! sequential vec on io_rank
189) !
190) ! Author: Glenn Hammond
191) ! Date: 10/12/12
192) !
193)
194) use Option_module
195) use Realization_Subsurface_class
196) use Grid_module
197) use Discretization_module
198)
199) implicit none
200)
201) #include "petsc/finclude/petscis.h"
202) #include "petsc/finclude/petscis.h90"
203) #include "petsc/finclude/petscviewer.h"
204)
205) type(regression_type), pointer :: regression
206) class(realization_subsurface_type) :: realization
207)
208) IS :: is_petsc
209) PetscInt, allocatable :: int_array(:)
210) PetscInt :: i, upper_bound, lower_bound, count, temp_int
211) PetscInt :: local_id
212) PetscReal, pointer :: vec_ptr(:)
213) Vec :: temp_vec
214) VecScatter :: temp_scatter
215) IS :: temp_is
216) PetscViewer :: viewer
217) PetscErrorCode :: ierr
218)
219) type(grid_type), pointer :: grid
220) type(option_type), pointer :: option
221)
222) if (.not.associated(regression)) return
223)
224) grid => realization%patch%grid
225) option => realization%option
226)
227) ! natural cell ids
228) if (associated(regression%natural_cell_ids)) then
229) ! ensure that natural ids are within problem domain
230) if (maxval(regression%natural_cell_ids) > grid%nmax) then
231) option%io_buffer = 'Natural IDs outside problem domain requested ' // &
232) 'for regression output. Removing non-existent IDs.'
233) call printWrnMsg(option)
234) count = 0
235) allocate(int_array(size(regression%natural_cell_ids)))
236) int_array = 0
237) do i = 1, size(regression%natural_cell_ids)
238) if (regression%natural_cell_ids(i) <= grid%nmax) then
239) count = count + 1
240) int_array(count) = regression%natural_cell_ids(i)
241) endif
242) enddo
243) ! reallocate array
244) deallocate(regression%natural_cell_ids)
245) allocate(regression%natural_cell_ids(count))
246) !geh: Since natural_cell_ids and int_array may now be of different sizes,
247) ! we need to be explicit about the values to copy. gfortran has
248) ! issues with this while Intel figures it out. Better to be explicit.
249) regression%natural_cell_ids = int_array(1:count)
250) deallocate(int_array)
251) endif
252) call VecCreate(PETSC_COMM_SELF,regression%natural_cell_id_vec, &
253) ierr);CHKERRQ(ierr)
254) if (option%myrank == option%io_rank) then
255) call VecSetSizes(regression%natural_cell_id_vec, &
256) size(regression%natural_cell_ids), &
257) PETSC_DECIDE,ierr);CHKERRQ(ierr)
258) else
259) call VecSetSizes(regression%natural_cell_id_vec,0, &
260) PETSC_DECIDE,ierr);CHKERRQ(ierr)
261) endif
262) call VecSetFromOptions(regression%natural_cell_id_vec,ierr);CHKERRQ(ierr)
263)
264) if (option%myrank == option%io_rank) then
265) count = size(regression%natural_cell_ids)
266) ! determine how many of the natural cell ids are local
267) allocate(int_array(count))
268) int_array = regression%natural_cell_ids
269) ! convert to zero based
270) int_array = int_array - 1
271) else
272) count = 0
273) allocate(int_array(count))
274) endif
275) call DiscretAOApplicationToPetsc(realization%discretization,int_array)
276)
277) ! create IS for global petsc cell ids
278) call ISCreateGeneral(option%mycomm,count,int_array,PETSC_COPY_VALUES, &
279) is_petsc,ierr);CHKERRQ(ierr)
280) deallocate(int_array)
281)
282) #ifdef REGRESSION_DEBUG
283) call PetscViewerASCIIOpen(option%mycomm, &
284) 'is_petsc_natural_cell_id.out', &
285) viewer,ierr);CHKERRQ(ierr)
286) call ISView(is_petsc,viewer,ierr);CHKERRQ(ierr)
287) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
288) #endif
289)
290) ! create scatter context
291) call VecScatterCreate(realization%field%work,is_petsc, &
292) regression%natural_cell_id_vec,PETSC_NULL_OBJECT, &
293) regression%scatter_natural_cell_id_gtos, &
294) ierr);CHKERRQ(ierr)
295)
296) call ISDestroy(is_petsc,ierr);CHKERRQ(ierr)
297)
298) #ifdef REGRESSION_DEBUG
299) call PetscViewerASCIIOpen(option%mycomm, &
300) 'regression_scatter_nat_cell_ids.out',viewer, &
301) ierr);CHKERRQ(ierr)
302) call VecScatterView(regression%scatter_natural_cell_id_gtos,viewer, &
303) ierr);CHKERRQ(ierr)
304) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
305) #endif
306)
307) endif
308)
309) if (regression%num_cells_per_process > 0) then
310) ! determine minimum number of cells per process
311) i = grid%nlmax
312) call MPI_Allreduce(i,count,ONE_INTEGER_MPI,MPIU_INTEGER,MPI_MIN, &
313) option%mycomm,ierr)
314) if (count < regression%num_cells_per_process) then
315) option%io_buffer = 'Number of cells per process for Regression ' // &
316) 'exceeds minimum number of cells per process. Truncating.'
317) call printMsg(option)
318) regression%num_cells_per_process = count
319) endif
320)
321) ! cells ids per processor
322) call VecCreate(PETSC_COMM_SELF,regression%cells_per_process_vec, &
323) ierr);CHKERRQ(ierr)
324) if (option%myrank == option%io_rank) then
325) call VecSetSizes(regression%cells_per_process_vec, &
326) regression%num_cells_per_process*option%mycommsize, &
327) PETSC_DECIDE,ierr);CHKERRQ(ierr)
328) else
329) call VecSetSizes(regression%cells_per_process_vec,ZERO_INTEGER, &
330) PETSC_DECIDE,ierr);CHKERRQ(ierr)
331) endif
332) call VecSetFromOptions(regression%cells_per_process_vec, &
333) ierr);CHKERRQ(ierr)
334)
335) ! create temporary vec to transfer down ids of cells
336) call VecCreate(option%mycomm,temp_vec,ierr);CHKERRQ(ierr)
337) call VecSetSizes(temp_vec, &
338) regression%num_cells_per_process, &
339) PETSC_DECIDE,ierr);CHKERRQ(ierr)
340) call VecSetFromOptions(temp_vec,ierr);CHKERRQ(ierr)
341)
342) ! calculate interval
343) call VecGetArrayF90(temp_vec,vec_ptr,ierr);CHKERRQ(ierr)
344) temp_int = grid%nlmax / regression%num_cells_per_process
345) do i = 1, regression%num_cells_per_process
346) vec_ptr(i) = temp_int*(i-1) + 1 + grid%global_offset
347) enddo
348) call VecRestoreArrayF90(temp_vec,vec_ptr,ierr);CHKERRQ(ierr)
349)
350) ! create temporary scatter to transfer values to io_rank
351) if (option%myrank == option%io_rank) then
352) count = option%mycommsize*regression%num_cells_per_process
353) ! determine how many of the natural cell ids are local
354) allocate(int_array(count))
355) do i = 1, count
356) int_array(i) = i
357) enddo
358) ! convert to zero based
359) int_array = int_array - 1
360) else
361) count = 0
362) allocate(int_array(count))
363) endif
364) call ISCreateGeneral(option%mycomm,count, &
365) int_array,PETSC_COPY_VALUES,temp_is, &
366) ierr);CHKERRQ(ierr)
367)
368) call VecScatterCreate(temp_vec,temp_is, &
369) regression%cells_per_process_vec,PETSC_NULL_OBJECT, &
370) temp_scatter,ierr);CHKERRQ(ierr)
371) call ISDestroy(temp_is,ierr);CHKERRQ(ierr)
372)
373) ! scatter ids to io_rank
374) call VecScatterBegin(temp_scatter,temp_vec, &
375) regression%cells_per_process_vec, &
376) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
377) call VecScatterEnd(temp_scatter,temp_vec, &
378) regression%cells_per_process_vec, &
379) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
380) call VecScatterDestroy(temp_scatter,ierr);CHKERRQ(ierr)
381) call VecDestroy(temp_vec,ierr);CHKERRQ(ierr)
382)
383) ! transfer cell ids into array for creating new scatter
384) if (option%myrank == option%io_rank) then
385) count = option%mycommsize*regression%num_cells_per_process
386) call VecGetArrayF90(regression%cells_per_process_vec,vec_ptr, &
387) ierr);CHKERRQ(ierr)
388) do i = 1, count
389) int_array(i) = int(vec_ptr(i)+0.1d0) ! tolerance to ensure int value
390) enddo
391) call VecRestoreArrayF90(regression%cells_per_process_vec,vec_ptr, &
392) ierr);CHKERRQ(ierr)
393) ! convert to zero based
394) int_array = int_array - 1
395) endif
396)
397) call ISCreateGeneral(option%mycomm,count, &
398) int_array,PETSC_COPY_VALUES,is_petsc, &
399) ierr);CHKERRQ(ierr)
400) deallocate(int_array)
401)
402) #ifdef REGRESSION_DEBUG
403) call PetscViewerASCIIOpen(option%mycomm, &
404) 'is_petsc_cells_per_process.out', &
405) viewer,ierr);CHKERRQ(ierr)
406) call ISView(is_petsc,viewer,ierr);CHKERRQ(ierr)
407) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
408) #endif
409)
410) call VecScatterCreate(realization%field%work,is_petsc, &
411) regression%cells_per_process_vec, &
412) PETSC_NULL_OBJECT, &
413) regression%scatter_cells_per_process_gtos, &
414) ierr);CHKERRQ(ierr)
415) call ISDestroy(is_petsc,ierr);CHKERRQ(ierr)
416)
417) #ifdef REGRESSION_DEBUG
418) call PetscViewerASCIIOpen(option%mycomm, &
419) 'regression_scatter_cells_per_process.out', &
420) viewer,ierr);CHKERRQ(ierr)
421) call VecScatterView(regression%scatter_cells_per_process_gtos,viewer, &
422) ierr);CHKERRQ(ierr)
423) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
424) #endif
425)
426) ! fill in natural ids of these cells on the io_rank
427) if (option%myrank == option%io_rank) then
428) allocate(regression%cells_per_process_natural_ids( &
429) regression%num_cells_per_process*option%mycommsize))
430) endif
431)
432) call VecGetArrayF90(realization%field%work,vec_ptr,ierr);CHKERRQ(ierr)
433) do local_id = 1, grid%nlmax
434) vec_ptr(local_id) = grid%nG2A(grid%nL2G(local_id))
435) enddo
436) call VecRestoreArrayF90(realization%field%work,vec_ptr,ierr);CHKERRQ(ierr)
437)
438) call VecScatterBegin(regression%scatter_cells_per_process_gtos, &
439) realization%field%work, &
440) regression%cells_per_process_vec, &
441) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
442) call VecScatterEnd(regression%scatter_cells_per_process_gtos, &
443) realization%field%work, &
444) regression%cells_per_process_vec, &
445) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
446)
447) if (option%myrank == option%io_rank) then
448) call VecGetArrayF90(regression%cells_per_process_vec,vec_ptr, &
449) ierr);CHKERRQ(ierr)
450) regression%cells_per_process_natural_ids(:) = int(vec_ptr(:)+0.1)
451) call VecRestoreArrayF90(regression%cells_per_process_vec,vec_ptr, &
452) ierr);CHKERRQ(ierr)
453) endif
454)
455) endif
456)
457) end subroutine RegressionCreateMapping
458)
459) ! ************************************************************************** !
460)
461) subroutine RegressionOutput(regression,realization,flow_timestepper, &
462) tran_timestepper)
463) !
464) ! Prints regression output through the io_rank
465) !
466) ! Author: Glenn Hammond
467) ! Date: 10/12/12
468) !
469)
470) use Realization_Subsurface_class
471) use Timestepper_BE_class
472) use Option_module
473) use Discretization_module
474) use Output_module
475) use Output_Aux_module
476) use Output_Common_module, only : OutputGetCellCenteredVelocities, &
477) OutputGetVarFromArray
478)
479) implicit none
480)
481) type(regression_type), pointer :: regression
482) class(realization_subsurface_type) :: realization
483) ! these must be pointers as they can be null
484) class(timestepper_BE_type), pointer :: flow_timestepper
485) class(timestepper_BE_type), pointer :: tran_timestepper
486)
487) character(len=MAXSTRINGLENGTH) :: string
488) Vec :: global_vec
489) Vec :: global_vec_vx,global_vec_vy,global_vec_vz
490) Vec :: x_vel_natural, y_vel_natural, z_vel_natural
491) Vec :: x_vel_process, y_vel_process, z_vel_process
492) PetscInt :: ivar, isubvar
493) type(option_type), pointer :: option
494) type(output_variable_type), pointer :: cur_variable
495) PetscReal, pointer :: vec_ptr(:), y_ptr(:), z_ptr(:)
496) PetscInt :: i
497) PetscInt :: iphase
498) PetscReal :: r_norm, x_norm
499) PetscReal :: max, min, mean
500) PetscErrorCode :: ierr
501)
502) if (.not.associated(regression)) return
503)
504) option => realization%option
505)
506) if (option%myrank == option%io_rank) then
507) string = trim(option%global_prefix) // &
508) trim(option%group_prefix) // &
509) '.regression'
510) option%io_buffer = '--> write regression output file: ' // trim(string)
511) call printMsg(option)
512) open(unit=OUTPUT_UNIT,file=string,action="write")
513) endif
514)
515) call DiscretizationCreateVector(realization%discretization,ONEDOF, &
516) global_vec,GLOBAL,option)
517) call DiscretizationDuplicateVector(realization%discretization,global_vec,global_vec_vx)
518) call DiscretizationDuplicateVector(realization%discretization,global_vec,global_vec_vy)
519) call DiscretizationDuplicateVector(realization%discretization,global_vec,global_vec_vz)
520)
521) cur_variable => realization%output_option%output_snap_variable_list%first
522) do
523) if (.not.associated(cur_variable)) exit
524)
525) ivar = cur_variable%ivar
526) isubvar = cur_variable%isubvar
527)
528) call OutputGetVarFromArray(realization,global_vec,ivar,isubvar)
529)
530) call VecMax(global_vec,PETSC_NULL_INTEGER,max,ierr);CHKERRQ(ierr)
531) call VecMin(global_vec,PETSC_NULL_INTEGER,min,ierr);CHKERRQ(ierr)
532) call VecSum(global_vec,mean,ierr);CHKERRQ(ierr)
533) mean = mean / realization%patch%grid%nmax
534)
535) ! list of natural ids
536) if (associated(regression%natural_cell_ids)) then
537) call VecScatterBegin(regression%scatter_natural_cell_id_gtos, &
538) global_vec, &
539) regression%natural_cell_id_vec, &
540) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
541) call VecScatterEnd(regression%scatter_natural_cell_id_gtos, &
542) global_vec, &
543) regression%natural_cell_id_vec, &
544) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
545) endif
546) if (regression%num_cells_per_process > 0) then
547) ! cells per process
548) call VecScatterBegin(regression%scatter_cells_per_process_gtos, &
549) global_vec, &
550) regression%cells_per_process_vec, &
551) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
552) call VecScatterEnd(regression%scatter_cells_per_process_gtos, &
553) global_vec, &
554) regression%cells_per_process_vec, &
555) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
556) endif
557)
558) 100 format(i9,': ',es21.13)
559) 101 format(i9,': ',i9)
560)
561) if (option%myrank == option%io_rank) then
562) string = OutputVariableToCategoryString(cur_variable%icategory)
563) write(OUTPUT_UNIT,'(''-- '',a,'': '',a,'' --'')') &
564) trim(string), trim(cur_variable%name)
565)
566) ! max, min, mean
567) if (cur_variable%iformat == 0) then
568) write(OUTPUT_UNIT,'(6x,''Max: '',es21.13)') max
569) write(OUTPUT_UNIT,'(6x,''Min: '',es21.13)') min
570) else
571) write(OUTPUT_UNIT,'(6x,''Max: '',i9)') int(max)
572) write(OUTPUT_UNIT,'(6x,''Min: '',i9)') int(min)
573) endif
574) write(OUTPUT_UNIT,'(5x,''Mean: '',es21.13)') mean
575)
576) ! natural cell ids
577) if (associated(regression%natural_cell_ids)) then
578) if (size(regression%natural_cell_ids) > 0) then
579) call VecGetArrayF90(regression%natural_cell_id_vec,vec_ptr, &
580) ierr);CHKERRQ(ierr)
581) if (cur_variable%iformat == 0) then
582) do i = 1, size(regression%natural_cell_ids)
583) write(OUTPUT_UNIT,100) &
584) regression%natural_cell_ids(i),vec_ptr(i)
585) enddo
586) else
587) do i = 1, size(regression%natural_cell_ids)
588) write(OUTPUT_UNIT,101) &
589) regression%natural_cell_ids(i),int(vec_ptr(i))
590) enddo
591) endif
592) call VecRestoreArrayF90(regression%natural_cell_id_vec,vec_ptr, &
593) ierr);CHKERRQ(ierr)
594) endif
595) endif
596)
597) ! cell ids per process
598) if (regression%num_cells_per_process > 0) then
599) call VecGetArrayF90(regression%cells_per_process_vec,vec_ptr, &
600) ierr);CHKERRQ(ierr)
601) if (cur_variable%iformat == 0) then
602) do i = 1, regression%num_cells_per_process*option%mycommsize
603) write(OUTPUT_UNIT,100) &
604) regression%cells_per_process_natural_ids(i),vec_ptr(i)
605) enddo
606) else
607) do i = 1, regression%num_cells_per_process*option%mycommsize
608) write(OUTPUT_UNIT,101) &
609) regression%cells_per_process_natural_ids(i),int(vec_ptr(i))
610) enddo
611) endif
612) call VecRestoreArrayF90(regression%cells_per_process_vec,vec_ptr, &
613) ierr);CHKERRQ(ierr)
614) endif
615) endif
616)
617) cur_variable => cur_variable%next
618) enddo
619)
620) ! velocities
621) if ((realization%output_option%print_tecplot_vel_cent .or. &
622) realization%output_option%print_hdf5_vel_cent) .and. &
623) option%nflowdof > 0) then
624) if (associated(regression%natural_cell_ids)) then
625) call VecDuplicate(regression%natural_cell_id_vec,x_vel_natural, &
626) ierr);CHKERRQ(ierr)
627) call VecDuplicate(regression%natural_cell_id_vec,y_vel_natural, &
628) ierr);CHKERRQ(ierr)
629) call VecDuplicate(regression%natural_cell_id_vec,z_vel_natural, &
630) ierr);CHKERRQ(ierr)
631) call VecZeroEntries(x_vel_natural,ierr);CHKERRQ(ierr)
632) call VecZeroEntries(y_vel_natural,ierr);CHKERRQ(ierr)
633) call VecZeroEntries(z_vel_natural,ierr);CHKERRQ(ierr)
634) endif
635) if (regression%num_cells_per_process > 0) then
636) call VecDuplicate(regression%cells_per_process_vec,x_vel_process, &
637) ierr);CHKERRQ(ierr)
638) call VecDuplicate(regression%cells_per_process_vec,y_vel_process, &
639) ierr);CHKERRQ(ierr)
640) call VecDuplicate(regression%cells_per_process_vec,z_vel_process, &
641) ierr);CHKERRQ(ierr)
642) call VecZeroEntries(x_vel_process,ierr);CHKERRQ(ierr)
643) call VecZeroEntries(y_vel_process,ierr);CHKERRQ(ierr)
644) call VecZeroEntries(z_vel_process,ierr);CHKERRQ(ierr)
645) endif
646)
647) do iphase = 1, option%nphase
648) if (associated(regression%natural_cell_ids) .or. &
649) regression%num_cells_per_process > 0) then
650)
651) if (iphase == 1) then
652) string = 'LIQUID'
653) else
654) string = 'GAS'
655) endif
656) if (option%myrank == option%io_rank) then
657) write(OUTPUT_UNIT,'(''-- GENERIC: '',a,'' VELOCITY ['',a, &
658) &''] --'')') &
659) trim(string), 'm/' // trim(realization%output_option%tunit)
660) endif
661)
662) ! X
663) call OutputGetCellCenteredVelocities(realization,global_vec_vx, &
664) global_vec_vy,global_vec_vz, &
665) iphase)
666) if (associated(regression%natural_cell_ids)) then
667) call VecScatterBegin(regression%scatter_natural_cell_id_gtos, &
668) global_vec_vx,x_vel_natural,INSERT_VALUES, &
669) SCATTER_FORWARD,ierr);CHKERRQ(ierr)
670) call VecScatterEnd(regression%scatter_natural_cell_id_gtos, &
671) global_vec_vx,x_vel_natural,INSERT_VALUES, &
672) SCATTER_FORWARD,ierr);CHKERRQ(ierr)
673) endif
674) if (regression%num_cells_per_process > 0) then
675) call VecScatterBegin(regression%scatter_cells_per_process_gtos, &
676) global_vec_vx,x_vel_process,INSERT_VALUES, &
677) SCATTER_FORWARD,ierr);CHKERRQ(ierr)
678) call VecScatterEnd(regression%scatter_cells_per_process_gtos, &
679) global_vec_vx,x_vel_process,INSERT_VALUES, &
680) SCATTER_FORWARD,ierr);CHKERRQ(ierr)
681) endif
682) ! Y
683) if (associated(regression%natural_cell_ids)) then
684) call VecScatterBegin(regression%scatter_natural_cell_id_gtos, &
685) global_vec_vy,y_vel_natural,INSERT_VALUES, &
686) SCATTER_FORWARD,ierr);CHKERRQ(ierr)
687) call VecScatterEnd(regression%scatter_natural_cell_id_gtos, &
688) global_vec_vy,y_vel_natural,INSERT_VALUES, &
689) SCATTER_FORWARD,ierr);CHKERRQ(ierr)
690) endif
691) if (regression%num_cells_per_process > 0) then
692) call VecScatterBegin(regression%scatter_cells_per_process_gtos, &
693) global_vec_vy,y_vel_process,INSERT_VALUES, &
694) SCATTER_FORWARD,ierr);CHKERRQ(ierr)
695) call VecScatterEnd(regression%scatter_cells_per_process_gtos, &
696) global_vec_vy,y_vel_process,INSERT_VALUES, &
697) SCATTER_FORWARD,ierr);CHKERRQ(ierr)
698) endif
699) ! Z
700) if (associated(regression%natural_cell_ids)) then
701) call VecScatterBegin(regression%scatter_natural_cell_id_gtos, &
702) global_vec_vz,z_vel_natural,INSERT_VALUES, &
703) SCATTER_FORWARD,ierr);CHKERRQ(ierr)
704) call VecScatterEnd(regression%scatter_natural_cell_id_gtos, &
705) global_vec_vz,z_vel_natural,INSERT_VALUES, &
706) SCATTER_FORWARD,ierr);CHKERRQ(ierr)
707) endif
708) if (regression%num_cells_per_process > 0) then
709) call VecScatterBegin(regression%scatter_cells_per_process_gtos, &
710) global_vec_vz,z_vel_process,INSERT_VALUES, &
711) SCATTER_FORWARD,ierr);CHKERRQ(ierr)
712) call VecScatterEnd(regression%scatter_cells_per_process_gtos, &
713) global_vec_vz,z_vel_process,INSERT_VALUES, &
714) SCATTER_FORWARD,ierr);CHKERRQ(ierr)
715) endif
716)
717) 104 format(i9,': ',3es21.13)
718)
719) ! natural cell ids
720) if (option%myrank == option%io_rank) then
721) if (associated(regression%natural_cell_ids)) then
722) if (size(regression%natural_cell_ids) > 0) then
723) call VecGetArrayF90(x_vel_natural,vec_ptr,ierr);CHKERRQ(ierr)
724) call VecGetArrayF90(y_vel_natural,y_ptr,ierr);CHKERRQ(ierr)
725) call VecGetArrayF90(z_vel_natural,z_ptr,ierr);CHKERRQ(ierr)
726) do i = 1, size(regression%natural_cell_ids)
727) write(OUTPUT_UNIT,104) &
728) regression%natural_cell_ids(i),vec_ptr(i),y_ptr(i),z_ptr(i)
729) enddo
730) call VecRestoreArrayF90(x_vel_natural,vec_ptr, &
731) ierr);CHKERRQ(ierr)
732) call VecRestoreArrayF90(y_vel_natural,y_ptr,ierr);CHKERRQ(ierr)
733) call VecRestoreArrayF90(z_vel_natural,z_ptr,ierr);CHKERRQ(ierr)
734) endif
735) endif
736)
737) ! cell ids per process
738) if (regression%num_cells_per_process > 0) then
739) call VecGetArrayF90(x_vel_process,vec_ptr,ierr);CHKERRQ(ierr)
740) call VecGetArrayF90(y_vel_process,y_ptr,ierr);CHKERRQ(ierr)
741) call VecGetArrayF90(z_vel_process,z_ptr,ierr);CHKERRQ(ierr)
742) do i = 1, regression%num_cells_per_process*option%mycommsize
743) write(OUTPUT_UNIT,104) &
744) regression%cells_per_process_natural_ids(i),vec_ptr(i), &
745) y_ptr(i),z_ptr(i)
746) enddo
747) call VecRestoreArrayF90(x_vel_process,vec_ptr,ierr);CHKERRQ(ierr)
748) call VecRestoreArrayF90(y_vel_process,y_ptr,ierr);CHKERRQ(ierr)
749) call VecRestoreArrayF90(z_vel_process,z_ptr,ierr);CHKERRQ(ierr)
750) endif
751) endif
752) endif
753) enddo
754)
755) if (associated(regression%natural_cell_ids)) then
756) call VecDestroy(x_vel_natural,ierr);CHKERRQ(ierr)
757) call VecDestroy(y_vel_natural,ierr);CHKERRQ(ierr)
758) call VecDestroy(z_vel_natural,ierr);CHKERRQ(ierr)
759) endif
760) if (regression%num_cells_per_process > 0) then
761) call VecDestroy(x_vel_process,ierr);CHKERRQ(ierr)
762) call VecDestroy(y_vel_process,ierr);CHKERRQ(ierr)
763) call VecDestroy(z_vel_process,ierr);CHKERRQ(ierr)
764) endif
765) endif ! option%nflowdof > 0
766)
767) call VecDestroy(global_vec,ierr);CHKERRQ(ierr)
768) call VecDestroy(global_vec_vx,ierr);CHKERRQ(ierr)
769) call VecDestroy(global_vec_vy,ierr);CHKERRQ(ierr)
770) call VecDestroy(global_vec_vz,ierr);CHKERRQ(ierr)
771)
772) 102 format(i12)
773) 103 format(es21.13)
774)
775) ! timestep, newton iteration, solver iteration output
776) if (associated(flow_timestepper)) then
777) call VecNorm(realization%field%flow_xx,NORM_2,x_norm,ierr);CHKERRQ(ierr)
778) call VecNorm(realization%field%flow_r,NORM_2,r_norm,ierr);CHKERRQ(ierr)
779) if (option%myrank == option%io_rank) then
780) write(OUTPUT_UNIT,'(''-- SOLUTION: Flow --'')')
781) write(OUTPUT_UNIT,'('' Time (seconds): '',es21.13)') &
782) flow_timestepper%cumulative_solver_time
783) write(OUTPUT_UNIT,'('' Time Steps: '',i12)') flow_timestepper%steps
784) write(OUTPUT_UNIT,'('' Newton Iterations: '',i12)') &
785) flow_timestepper%cumulative_newton_iterations
786) write(OUTPUT_UNIT,'('' Solver Iterations: '',i12)') &
787) flow_timestepper%cumulative_linear_iterations
788) write(OUTPUT_UNIT,'('' Time Step Cuts: '',i12)') &
789) flow_timestepper%cumulative_time_step_cuts
790) write(OUTPUT_UNIT,'('' Solution 2-Norm: '',es21.13)') x_norm
791) write(OUTPUT_UNIT,'('' Residual 2-Norm: '',es21.13)') r_norm
792) endif
793) endif
794) if (associated(tran_timestepper)) then
795) call VecNorm(realization%field%tran_xx,NORM_2,x_norm,ierr);CHKERRQ(ierr)
796) if (option%transport%reactive_transport_coupling == GLOBAL_IMPLICIT) then
797) call VecNorm(realization%field%tran_r,NORM_2,r_norm,ierr);CHKERRQ(ierr)
798) endif
799) if (option%myrank == option%io_rank) then
800) write(OUTPUT_UNIT,'(''-- SOLUTION: Transport --'')')
801) write(OUTPUT_UNIT,'('' Time (seconds): '',es21.13)') &
802) tran_timestepper%cumulative_solver_time
803) write(OUTPUT_UNIT,'('' Time Steps: '',i12)') tran_timestepper%steps
804) write(OUTPUT_UNIT,'('' Newton Iterations: '',i12)') &
805) tran_timestepper%cumulative_newton_iterations
806) write(OUTPUT_UNIT,'('' Solver Iterations: '',i12)') &
807) tran_timestepper%cumulative_linear_iterations
808) write(OUTPUT_UNIT,'('' Time Step Cuts: '',i12)') &
809) tran_timestepper%cumulative_time_step_cuts
810) write(OUTPUT_UNIT,'('' Solution 2-Norm: '',es21.13)') x_norm
811) if (option%transport%reactive_transport_coupling == GLOBAL_IMPLICIT) then
812) write(OUTPUT_UNIT,'('' Residual 2-Norm: '',es21.13)') r_norm
813) endif
814) endif
815) endif
816)
817) close(OUTPUT_UNIT)
818)
819) end subroutine RegressionOutput
820)
821) ! ************************************************************************** !
822)
823) recursive subroutine RegressionVariableDestroy(regression_variable)
824) !
825) ! Destroys a regression variable object
826) !
827) ! Author: Glenn Hammond
828) ! Date: 10/11/12
829) !
830)
831) implicit none
832)
833) type(regression_variable_type), pointer :: regression_variable
834)
835) if (.not.associated(regression_variable)) return
836)
837) call RegressionVariableDestroy(regression_variable%next)
838)
839) deallocate(regression_variable)
840) nullify(regression_variable)
841)
842) end subroutine RegressionVariableDestroy
843)
844) ! ************************************************************************** !
845)
846) subroutine RegressionDestroy(regression)
847) !
848) ! Destroys a regression object
849) !
850) ! Author: Glenn Hammond
851) ! Date: 10/11/12
852) !
853)
854) use Utility_module
855)
856) implicit none
857)
858) type(regression_type), pointer :: regression
859)
860) PetscErrorCode :: ierr
861)
862) if (.not.associated(regression)) return
863)
864) call RegressionVariableDestroy(regression%variable_list)
865) call DeallocateArray(regression%natural_cell_ids)
866) regression%num_cells_per_process = 0
867) call DeallocateArray(regression%cells_per_process_natural_ids)
868) if (regression%natural_cell_id_vec /= 0) then
869) call VecDestroy(regression%natural_cell_id_vec,ierr);CHKERRQ(ierr)
870) endif
871) if (regression%cells_per_process_vec /= 0) then
872) call VecDestroy(regression%cells_per_process_vec,ierr);CHKERRQ(ierr)
873) endif
874) if (regression%scatter_natural_cell_id_gtos /= 0) then
875) call VecScatterDestroy(regression%scatter_natural_cell_id_gtos, &
876) ierr);CHKERRQ(ierr)
877) endif
878) if (regression%scatter_cells_per_process_gtos /= 0) then
879) call VecScatterDestroy(regression%scatter_cells_per_process_gtos, &
880) ierr);CHKERRQ(ierr)
881) endif
882)
883) deallocate(regression)
884) nullify(regression)
885)
886) end subroutine RegressionDestroy
887)
888) end module Regression_module