timestepper_BE.F90 coverage: 94.44 %func 87.13 %block
1) module Timestepper_BE_class
2)
3) use Solver_module
4) use Convergence_module
5) use Timestepper_Base_class
6)
7) use PFLOTRAN_Constants_module
8)
9) implicit none
10)
11) private
12)
13) #include "petsc/finclude/petscsys.h"
14)
15) type, public, extends(timestepper_base_type) :: timestepper_BE_type
16)
17) PetscInt :: num_newton_iterations ! number of Newton iterations in a time step
18) PetscInt :: num_linear_iterations ! number of linear solver iterations in a time step
19) PetscInt :: cumulative_newton_iterations ! Total number of Newton iterations
20) PetscInt :: cumulative_linear_iterations ! Total number of linear iterations
21) PetscInt :: cumulative_wasted_linear_iterations
22)
23) PetscInt :: iaccel ! Accelerator index
24) ! An array of multiplicative factors that specify how to increase time step.
25) PetscReal, pointer :: tfac(:)
26) PetscInt :: ntfac ! size of tfac
27)
28) type(solver_type), pointer :: solver
29) type(convergence_context_type), pointer :: convergence_context
30)
31) contains
32)
33) procedure, public :: ReadInput => TimestepperBERead
34) procedure, public :: Init => TimestepperBEInit
35) ! procedure, public :: SetTargetTime => TimestepperBaseSetTargetTime
36) procedure, public :: StepDT => TimestepperBEStepDT
37) procedure, public :: UpdateDT => TimestepperBEUpdateDT
38) procedure, public :: CheckpointBinary => TimestepperBECheckpointBinary
39) procedure, public :: RestartBinary => TimestepperBERestartBinary
40) #if defined(PETSC_HAVE_HDF5)
41) procedure, public :: CheckpointHDF5 => TimestepperBECheckpointHDF5
42) procedure, public :: RestartHDF5 => TimestepperBERestartHDF5
43) #endif
44) procedure, public :: Reset => TimestepperBEReset
45) procedure, public :: PrintInfo => TimestepperBEPrintInfo
46) procedure, public :: InputRecord => TimestepperBEInputRecord
47) procedure, public :: FinalizeRun => TimestepperBEFinalizeRun
48) procedure, public :: Strip => TimestepperBEStrip
49) procedure, public :: Destroy => TimestepperBEDestroy
50)
51) end type timestepper_BE_type
52)
53) ! For checkpointing
54) type, public, extends(stepper_base_header_type) :: stepper_BE_header_type
55) PetscInt :: cumulative_newton_iterations
56) PetscInt :: cumulative_linear_iterations
57) PetscInt :: num_newton_iterations
58) end type stepper_BE_header_type
59)
60) interface PetscBagGetData
61) subroutine PetscBagGetData(bag,header,ierr)
62) import :: stepper_BE_header_type
63) implicit none
64) #include "petsc/finclude/petscbag.h"
65) PetscBag :: bag
66) class(stepper_BE_header_type), pointer :: header
67) PetscErrorCode :: ierr
68) end subroutine
69) end interface PetscBagGetData
70)
71) public :: TimestepperBECreate, TimestepperBEPrintInfo, &
72) TimestepperBEInit
73)
74) contains
75)
76) ! ************************************************************************** !
77)
78) function TimestepperBECreate()
79) !
80) ! Allocates and initializes a new Timestepper object
81) !
82) ! Author: Glenn Hammond
83) ! Date: 07/22/13
84) !
85)
86) implicit none
87)
88) class(timestepper_BE_type), pointer :: TimestepperBECreate
89)
90) class(timestepper_BE_type), pointer :: stepper
91)
92) allocate(stepper)
93) call stepper%Init()
94)
95) stepper%solver => SolverCreate()
96)
97) TimestepperBECreate => stepper
98)
99) end function TimestepperBECreate
100)
101) ! ************************************************************************** !
102)
103) subroutine TimestepperBEInit(this)
104) !
105) ! Allocates and initializes a new Timestepper object
106) !
107) ! Author: Glenn Hammond
108) ! Date: 07/22/13
109) !
110)
111) implicit none
112)
113) class(timestepper_BE_type) :: this
114)
115) call TimestepperBaseInit(this)
116)
117) this%num_newton_iterations = 0
118) this%num_linear_iterations = 0
119)
120) this%cumulative_newton_iterations = 0
121) this%cumulative_linear_iterations = 0
122) this%cumulative_wasted_linear_iterations = 0
123)
124) this%iaccel = 5
125) this%ntfac = 13
126) allocate(this%tfac(13))
127) this%tfac(1) = 2.0d0; this%tfac(2) = 2.0d0
128) this%tfac(3) = 2.0d0; this%tfac(4) = 2.0d0
129) this%tfac(5) = 2.0d0; this%tfac(6) = 1.8d0
130) this%tfac(7) = 1.6d0; this%tfac(8) = 1.4d0
131) this%tfac(9) = 1.2d0; this%tfac(10) = 1.0d0
132) this%tfac(11) = 1.0d0; this%tfac(12) = 1.0d0
133) this%tfac(13) = 1.0d0
134)
135) nullify(this%solver)
136) nullify(this%convergence_context)
137)
138) end subroutine TimestepperBEInit
139)
140) ! ************************************************************************** !
141)
142) subroutine TimestepperBERead(this,input,option)
143) !
144) ! Reads parameters associated with time stepper
145) !
146) ! Author: Glenn Hammond
147) ! Date: 07/22/13
148) !
149)
150) use Option_module
151) use String_module
152) use Input_Aux_module
153) use Utility_module
154)
155) implicit none
156)
157) class(timestepper_BE_type) :: this
158) type(input_type), pointer :: input
159) type(option_type) :: option
160)
161) character(len=MAXWORDLENGTH) :: keyword
162) character(len=MAXSTRINGLENGTH) :: string
163)
164) input%ierr = 0
165) do
166)
167) call InputReadPflotranString(input,option)
168)
169) if (InputCheckExit(input,option)) exit
170)
171) call InputReadWord(input,option,keyword,PETSC_TRUE)
172) call InputErrorMsg(input,option,'keyword','TIMESTEPPER_BE')
173) call StringToUpper(keyword)
174)
175) select case(trim(keyword))
176)
177) case('TS_ACCELERATION')
178) call InputReadInt(input,option,this%iaccel)
179) call InputDefaultMsg(input,option,'iaccel')
180)
181) case('DT_FACTOR')
182) string='time_step_factor'
183) call UtilityReadArray(this%tfac,NEG_ONE_INTEGER,string,input, &
184) option)
185) this%ntfac = size(this%tfac)
186)
187) case default
188) call TimestepperBaseProcessKeyword(this,input,option,keyword)
189) end select
190)
191) enddo
192)
193) this%solver%print_ekg = this%print_ekg
194)
195) end subroutine TimestepperBERead
196)
197) ! ************************************************************************** !
198)
199) subroutine TimestepperBEUpdateDT(this,process_model)
200) !
201) ! Updates time step
202) !
203) ! Author: Glenn Hammond
204) ! Date: 07/22/13
205) !
206)
207) use PM_Base_class
208)
209) implicit none
210)
211) class(timestepper_BE_type) :: this
212) class(pm_base_type) :: process_model
213)
214) PetscBool :: update_time_step
215)
216) update_time_step = PETSC_TRUE
217)
218) if (this%time_step_cut_flag) then
219) this%num_constant_time_steps = 1
220) else if (this%num_constant_time_steps > 0) then
221) ! otherwise, only increment if the constant time step counter was
222) ! initialized to 1
223) this%num_constant_time_steps = &
224) this%num_constant_time_steps + 1
225) endif
226)
227) ! num_constant_time_steps = 0: normal time stepping with growing steps
228) ! num_constant_time_steps > 0: restriction of constant time steps until
229) ! constant_time_step_threshold is met
230) if (this%num_constant_time_steps > &
231) this%constant_time_step_threshold) then
232) this%num_constant_time_steps = 0
233) else if (this%num_constant_time_steps > 0) then
234) ! do not increase time step size
235) update_time_step = PETSC_FALSE
236) endif
237)
238) if (update_time_step .and. this%iaccel /= 0) then
239)
240) call process_model%UpdateTimestep(this%dt, &
241) this%dt_min, &
242) this%dt_max, &
243) this%iaccel, &
244) this%num_newton_iterations, &
245) this%tfac)
246)
247) endif
248)
249) end subroutine TimestepperBEUpdateDT
250)
251) ! ************************************************************************** !
252)
253) subroutine TimestepperBEStepDT(this,process_model,stop_flag)
254) !
255) ! Steps forward one step in time
256) !
257) ! Author: Glenn Hammond
258) ! Date: 07/22/13
259) !
260)
261) use PM_Base_class
262) use Option_module
263) use Output_module, only : Output
264) use Output_EKG_module, only : IUNIT_EKG
265)
266) implicit none
267)
268) #include "petsc/finclude/petscvec.h"
269) #include "petsc/finclude/petscvec.h90"
270) #include "petsc/finclude/petscsnes.h"
271)
272) class(timestepper_BE_type) :: this
273) class(pm_base_type) :: process_model
274) PetscInt :: stop_flag
275)
276) SNESConvergedReason :: snes_reason
277) PetscInt :: icut
278)
279) type(solver_type), pointer :: solver
280) type(option_type), pointer :: option
281)
282) PetscLogDouble :: log_start_time
283) PetscLogDouble :: log_end_time
284) PetscInt :: num_newton_iterations
285) PetscInt :: num_linear_iterations
286) PetscInt :: num_linear_iterations2
287) PetscInt :: sum_newton_iterations
288) PetscInt :: sum_linear_iterations
289) PetscInt :: sum_wasted_linear_iterations
290) character(len=MAXWORDLENGTH) :: tunit
291) PetscReal :: tconv
292) PetscReal :: fnorm, inorm, scaled_fnorm
293) PetscBool :: snapshot_plot_flag, observation_plot_flag, massbal_plot_flag
294) Vec :: residual_vec
295) PetscErrorCode :: ierr
296)
297) solver => this%solver
298) option => process_model%option
299)
300) !geh: for debugging
301) #ifdef DEBUG
302) write(process_model%option%io_buffer,'(es12.5)') this%dt
303) process_model%option%io_buffer = 'StepperStepDT(' // &
304) trim(adjustl(process_model%option%io_buffer)) // ')'
305) call printMsg(process_model%option)
306) #endif
307)
308) tconv = process_model%output_option%tconv
309) tunit = process_model%output_option%tunit
310) sum_linear_iterations = 0
311) sum_wasted_linear_iterations = 0
312) sum_newton_iterations = 0
313) icut = 0
314)
315) option%dt = this%dt
316) option%time = this%target_time-this%dt
317)
318) call process_model%InitializeTimestep()
319)
320) do
321)
322) call process_model%PreSolve()
323)
324) call PetscTime(log_start_time, ierr);CHKERRQ(ierr)
325)
326) call SNESSolve(solver%snes,PETSC_NULL_OBJECT, &
327) process_model%solution_vec,ierr);CHKERRQ(ierr)
328)
329) call PetscTime(log_end_time, ierr);CHKERRQ(ierr)
330)
331) this%cumulative_solver_time = &
332) this%cumulative_solver_time + &
333) (log_end_time - log_start_time)
334)
335) call SNESGetIterationNumber(solver%snes,num_newton_iterations, &
336) ierr);CHKERRQ(ierr)
337) call SNESGetLinearSolveIterations(solver%snes,num_linear_iterations, &
338) ierr);CHKERRQ(ierr)
339) call SNESGetConvergedReason(solver%snes,snes_reason,ierr);CHKERRQ(ierr)
340)
341) sum_newton_iterations = sum_newton_iterations + num_newton_iterations
342) sum_linear_iterations = sum_linear_iterations + num_linear_iterations
343)
344) if (snes_reason <= 0 .or. .not. process_model%AcceptSolution()) then
345) sum_wasted_linear_iterations = sum_wasted_linear_iterations + &
346) num_linear_iterations
347) ! The Newton solver diverged, so try reducing the time step.
348) icut = icut + 1
349) this%time_step_cut_flag = PETSC_TRUE
350) ! if a cut occurs on the last time step, the stop_flag will have been
351) ! set to TS_STOP_END_SIMULATION. Set back to TS_CONTINUE to prevent
352) ! premature ending of simulation.
353) if (stop_flag /= TS_STOP_MAX_TIME_STEP) stop_flag = TS_CONTINUE
354)
355) if (icut > this%max_time_step_cuts .or. this%dt < this%dt_min) then
356)
357) write(option%io_buffer,'(" Stopping: Time step cut criteria &
358) &exceeded!")')
359) call printMsg(option)
360) write(option%io_buffer,'(" icut =",i3,", max_time_step_cuts=",i3)') &
361) icut,this%max_time_step_cuts
362) call printMsg(option)
363) write(option%io_buffer,'(" dt =",es15.7,", dt_min=",es15.7)') &
364) this%dt/tconv,this%dt_min/tconv
365) call printMsg(option)
366)
367) process_model%output_option%plot_name = 'flow_cut_to_failure'
368) snapshot_plot_flag = PETSC_TRUE
369) observation_plot_flag = PETSC_FALSE
370) massbal_plot_flag = PETSC_FALSE
371) call Output(process_model%realization_base,snapshot_plot_flag, &
372) observation_plot_flag,massbal_plot_flag)
373) stop_flag = TS_STOP_FAILURE
374) return
375) endif
376)
377) this%target_time = this%target_time - this%dt
378)
379) this%dt = 0.5d0 * this%dt
380)
381) write(option%io_buffer,'('' -> Cut time step: snes='',i3, &
382) & '' icut= '',i2,''['',i3,'']'','' t= '',1pe12.5, '' dt= '', &
383) & 1pe12.5)') snes_reason,icut,this%cumulative_time_step_cuts, &
384) option%time/tconv, &
385) this%dt/tconv
386) call printMsg(option)
387) if (snes_reason == SNES_DIVERGED_LINEAR_SOLVE) then
388) call KSPGetIterationNumber(solver%ksp,num_linear_iterations2, &
389) ierr);CHKERRQ(ierr)
390) sum_wasted_linear_iterations = sum_wasted_linear_iterations + &
391) num_linear_iterations2
392) sum_linear_iterations = sum_linear_iterations + num_linear_iterations2
393) call SolverLinearPrintFailedReason(solver,option)
394) endif
395)
396) this%target_time = this%target_time + this%dt
397) option%dt = this%dt
398) call process_model%TimeCut()
399)
400) else
401) ! The Newton solver converged, so we can exit.
402) exit
403) endif
404) enddo
405)
406) this%steps = this%steps + 1
407) this%cumulative_newton_iterations = &
408) this%cumulative_newton_iterations + sum_newton_iterations
409) this%cumulative_linear_iterations = &
410) this%cumulative_linear_iterations + sum_linear_iterations
411) this%cumulative_wasted_linear_iterations = &
412) this%cumulative_wasted_linear_iterations + sum_wasted_linear_iterations
413) this%cumulative_time_step_cuts = &
414) this%cumulative_time_step_cuts + icut
415)
416) this%num_newton_iterations = num_newton_iterations
417) this%num_linear_iterations = num_linear_iterations
418)
419) ! print screen output
420) call SNESGetFunction(solver%snes,residual_vec,PETSC_NULL_OBJECT, &
421) PETSC_NULL_INTEGER,ierr);CHKERRQ(ierr)
422) call VecNorm(residual_vec,NORM_2,fnorm,ierr);CHKERRQ(ierr)
423) call VecNorm(residual_vec,NORM_INFINITY,inorm,ierr);CHKERRQ(ierr)
424) if (option%print_screen_flag) then
425) write(*, '(/," Step ",i6," Time= ",1pe12.5," Dt= ",1pe12.5, &
426) & " [",a,"]", " snes_conv_reason: ",i4,/," newton = ",i3, &
427) & " [",i8,"]", " linear = ",i5," [",i10,"]"," cuts = ",i2, &
428) & " [",i4,"]")') &
429) this%steps, &
430) this%target_time/tconv, &
431) this%dt/tconv, &
432) trim(tunit),snes_reason,sum_newton_iterations, &
433) this%cumulative_newton_iterations,sum_linear_iterations, &
434) this%cumulative_linear_iterations,icut, &
435) this%cumulative_time_step_cuts
436)
437)
438) if (associated(process_model%realization_base%discretization%grid)) then
439) scaled_fnorm = fnorm/process_model%realization_base% &
440) discretization%grid%nmax
441) else
442) scaled_fnorm = fnorm
443) endif
444)
445) print *,' --> SNES Linear/Non-Linear Iterations = ', &
446) num_linear_iterations,' / ',num_newton_iterations
447) write(*,'(" --> SNES Residual: ",1p3e14.6)') fnorm, scaled_fnorm, inorm
448) endif
449) if (option%print_file_flag) then
450) write(option%fid_out, '(" Step ",i6," Time= ",1pe12.5," Dt= ",1pe12.5, &
451) & " [",a,"]"," snes_conv_reason: ",i4,/," newton = ",i3, &
452) & " [",i8,"]", " linear = ",i5," [",i10,"]"," cuts = ",i2," [",i4,"]")') &
453) this%steps, &
454) this%target_time/tconv, &
455) this%dt/tconv, &
456) trim(tunit),snes_reason,sum_newton_iterations, &
457) this%cumulative_newton_iterations,sum_linear_iterations, &
458) this%cumulative_linear_iterations,icut, &
459) this%cumulative_time_step_cuts
460) endif
461)
462) option%time = this%target_time
463) call process_model%FinalizeTimestep()
464)
465) if (this%print_ekg .and. OptionPrintToFile(option)) then
466) 100 format(a32," TIMESTEP ",i10,2es16.8,a,i3,i5,i3,i5,i5,i10)
467) write(IUNIT_EKG,100) trim(this%name), this%steps, this%target_time/tconv, &
468) this%dt/tconv, trim(tunit), &
469) icut, this%cumulative_time_step_cuts, &
470) sum_newton_iterations, this%cumulative_newton_iterations, &
471) sum_linear_iterations, this%cumulative_linear_iterations
472) endif
473)
474) if (option%print_screen_flag) print *, ""
475)
476) end subroutine TimestepperBEStepDT
477)
478) ! ************************************************************************** !
479)
480) subroutine TimestepperBECheckpointBinary(this,viewer,option)
481) !
482) ! Checkpoints parameters/variables associated with
483) ! a time stepper.
484) !
485) ! Author: Glenn Hammond
486) ! Date: 07/25/13
487) !
488)
489) use Option_module
490)
491) implicit none
492)
493) #include "petsc/finclude/petscviewer.h"
494) #include "petsc/finclude/petscbag.h"
495)
496) class(timestepper_BE_type) :: this
497) PetscViewer :: viewer
498) type(option_type) :: option
499)
500) class(stepper_BE_header_type), pointer :: header
501) type(stepper_BE_header_type) :: dummy_header
502) character(len=1),pointer :: dummy_char(:)
503) PetscBag :: bag
504) PetscSizeT :: bagsize
505) PetscErrorCode :: ierr
506)
507) bagsize = size(transfer(dummy_header,dummy_char))
508)
509) call PetscBagCreate(option%mycomm,bagsize,bag,ierr);CHKERRQ(ierr)
510) call PetscBagGetData(bag,header,ierr);CHKERRQ(ierr)
511) call TimestepperBERegisterHeader(this,bag,header)
512) call TimestepperBESetHeader(this,bag,header)
513) call PetscBagView(bag,viewer,ierr);CHKERRQ(ierr)
514) call PetscBagDestroy(bag,ierr);CHKERRQ(ierr)
515)
516) end subroutine TimestepperBECheckpointBinary
517)
518) ! ************************************************************************** !
519)
520) subroutine TimestepperBERegisterHeader(this,bag,header)
521) !
522) ! Register header entries.
523) !
524) ! Author: Glenn Hammond
525) ! Date: 07/30/13
526) !
527)
528) use Option_module
529)
530) implicit none
531)
532) #include "petsc/finclude/petscviewer.h"
533) #include "petsc/finclude/petscbag.h"
534)
535) class(timestepper_BE_type) :: this
536) class(stepper_BE_header_type) :: header
537) PetscBag :: bag
538)
539) PetscErrorCode :: ierr
540)
541) call PetscBagRegisterInt(bag,header%cumulative_newton_iterations,0, &
542) "cumulative_newton_iterations","", &
543) ierr);CHKERRQ(ierr)
544) call PetscBagRegisterInt(bag,header%cumulative_linear_iterations,0, &
545) "cumulative_linear_iterations","", &
546) ierr);CHKERRQ(ierr)
547) ! need to add cumulative wasted linear iterations
548) call PetscBagRegisterInt(bag,header%num_newton_iterations,0, &
549) "num_newton_iterations","",ierr);CHKERRQ(ierr)
550)
551) call TimestepperBaseRegisterHeader(this,bag,header)
552)
553) end subroutine TimestepperBERegisterHeader
554)
555) ! ************************************************************************** !
556)
557) subroutine TimestepperBESetHeader(this,bag,header)
558) !
559) ! Sets values in checkpoint header.
560) !
561) ! Author: Glenn Hammond
562) ! Date: 07/25/13
563) !
564)
565) use Option_module
566)
567) implicit none
568)
569) #include "petsc/finclude/petscviewer.h"
570) #include "petsc/finclude/petscbag.h"
571)
572) class(timestepper_BE_type) :: this
573) class(stepper_BE_header_type) :: header
574) PetscBag :: bag
575)
576) PetscErrorCode :: ierr
577)
578) header%cumulative_newton_iterations = this%cumulative_newton_iterations
579) header%cumulative_linear_iterations = this%cumulative_linear_iterations
580) header%num_newton_iterations = this%num_newton_iterations
581)
582) call TimestepperBaseSetHeader(this,bag,header)
583)
584) end subroutine TimestepperBESetHeader
585)
586) ! ************************************************************************** !
587)
588) subroutine TimestepperBERestartBinary(this,viewer,option)
589) !
590) ! Checkpoints parameters/variables associated with
591) ! a time stepper.
592) !
593) ! Author: Glenn Hammond
594) ! Date: 07/25/13
595) !
596)
597) use Option_module
598)
599) implicit none
600)
601) #include "petsc/finclude/petscviewer.h"
602) #include "petsc/finclude/petscbag.h"
603)
604) class(timestepper_BE_type) :: this
605) PetscViewer :: viewer
606) type(option_type) :: option
607)
608) class(stepper_BE_header_type), pointer :: header
609) type(stepper_BE_header_type) :: dummy_header
610) character(len=1),pointer :: dummy_char(:)
611) PetscBag :: bag
612) PetscSizeT :: bagsize
613) PetscErrorCode :: ierr
614)
615) bagsize = size(transfer(dummy_header,dummy_char))
616)
617) call PetscBagCreate(option%mycomm,bagsize,bag,ierr);CHKERRQ(ierr)
618) call PetscBagGetData(bag,header,ierr);CHKERRQ(ierr)
619) call TimestepperBERegisterHeader(this,bag,header)
620) call PetscBagLoad(viewer,bag,ierr);CHKERRQ(ierr)
621) call TimestepperBEGetHeader(this,header)
622) call PetscBagDestroy(bag,ierr);CHKERRQ(ierr)
623)
624) end subroutine TimestepperBERestartBinary
625)
626) ! ************************************************************************** !
627)
628) #if defined(PETSC_HAVE_HDF5)
629) subroutine TimestepperBECheckpointHDF5(this, chk_grp_id, option)
630) !
631) ! Checkpoints parameters/variables associated with
632) ! a time stepper.
633) !
634) ! Author: Gautam Bisht
635) ! Date: 07/30/15
636) !
637) use Option_module
638) use hdf5
639) use Checkpoint_module, only : CheckPointWriteIntDatasetHDF5
640) use Checkpoint_module, only : CheckPointWriteRealDatasetHDF5
641)
642) implicit none
643)
644) class(timestepper_BE_type) :: this
645) PetscInt :: chk_grp_id
646) type(option_type) :: option
647)
648) #if defined(SCORPIO_WRITE)
649) integer :: h5_chk_grp_id
650) integer, pointer :: dims(:)
651) integer, pointer :: start(:)
652) integer, pointer :: stride(:)
653) integer, pointer :: length(:)
654) integer :: timestepper_grp_id
655) #else
656) integer(HSIZE_T), pointer :: dims(:)
657) integer(HSIZE_T), pointer :: start(:)
658) integer(HSIZE_T), pointer :: stride(:)
659) integer(HSIZE_T), pointer :: length(:)
660) integer(HID_T) :: timestepper_grp_id
661) integer(HID_T) :: h5_chk_grp_id
662) #endif
663)
664) PetscMPIInt :: dataset_rank
665) character(len=MAXSTRINGLENGTH) :: dataset_name
666) character(len=MAXSTRINGLENGTH) :: string
667) PetscInt, pointer :: int_array(:)
668) PetscReal, pointer :: real_array(:)
669) PetscMPIInt :: hdf5_err
670)
671) string = "Timestepper"
672) call h5gcreate_f(chk_grp_id, string, timestepper_grp_id, hdf5_err, OBJECT_NAMELEN_DEFAULT_F)
673)
674) allocate(start(1))
675) allocate(dims(1))
676) allocate(length(1))
677) allocate(stride(1))
678) allocate(int_array(1))
679) allocate(real_array(1))
680)
681) dataset_rank = 1
682) dims(1) = ONE_INTEGER
683) start(1) = 0
684) length(1) = ONE_INTEGER
685) stride(1) = ONE_INTEGER
686)
687) dataset_name = "Cumulative_newton_iterations" // CHAR(0)
688) int_array(1) = this%cumulative_newton_iterations
689) call CheckPointWriteIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
690) dims, start, length, stride, int_array, option)
691)
692) dataset_name = "Cumulative_linear_iterations" // CHAR(0)
693) int_array(1) = this%cumulative_linear_iterations
694) call CheckPointWriteIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
695) dims, start, length, stride, int_array, option)
696)
697) dataset_name = "Num_newton_iterations" // CHAR(0)
698) int_array(1) = this%num_newton_iterations
699) call CheckPointWriteIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
700) dims, start, length, stride, int_array, option)
701)
702) dataset_name = "Time" // CHAR(0)
703) real_array(1) = this%target_time
704) call CheckPointWriteRealDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
705) dims, start, length, stride, real_array, option)
706)
707) dataset_name = "Dt" // CHAR(0)
708) real_array(1) = this%dt
709) call CheckPointWriteRealDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
710) dims, start, length, stride, real_array, option)
711)
712) dataset_name = "Prev_dt" // CHAR(0)
713) real_array(1) = this%prev_dt
714) call CheckPointWriteRealDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
715) dims, start, length, stride, real_array, option)
716)
717) dataset_name = "Num_steps" // CHAR(0)
718) int_array(1) = this%steps
719) call CheckPointWriteIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
720) dims, start, length, stride, int_array, option)
721)
722) dataset_name = "Cumulative_time_step_cuts" // CHAR(0)
723) int_array(1) = this%cumulative_time_step_cuts
724) call CheckPointWriteIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
725) dims, start, length, stride, int_array, option)
726)
727) dataset_name = "Num_constant_time_steps" // CHAR(0)
728) int_array(1) = this%num_constant_time_steps
729) call CheckPointWriteIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
730) dims, start, length, stride, int_array, option)
731)
732) dataset_name = "Num_contig_revert_due_to_sync" // CHAR(0)
733) int_array(1) = this%num_contig_revert_due_to_sync
734) call CheckPointWriteIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
735) dims, start, length, stride, int_array, option)
736)
737) dataset_name = "Revert_dt" // CHAR(0)
738) int_array(1) = ZERO_INTEGER
739) if (this%revert_dt) int_array(1) = ONE_INTEGER
740) call CheckPointWriteIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
741) dims, start, length, stride, int_array, option)
742)
743) call h5gclose_f(timestepper_grp_id, hdf5_err)
744)
745) deallocate(start)
746) deallocate(dims)
747) deallocate(length)
748) deallocate(stride)
749) deallocate(int_array)
750) deallocate(real_array)
751)
752) end subroutine TimestepperBECheckpointHDF5
753)
754) ! ************************************************************************** !
755)
756) subroutine TimestepperBERestartHDF5(this, chk_grp_id, option)
757) !
758) ! Restarts parameters/variables associated with
759) ! a time stepper.
760) !
761) ! Author: Gautam Bisht
762) ! Date: 08/16/15
763) !
764) use Option_module
765) use hdf5
766) use Checkpoint_module, only : CheckPointReadIntDatasetHDF5
767) use Checkpoint_module, only : CheckPointReadRealDatasetHDF5
768)
769) implicit none
770)
771) class(timestepper_BE_type) :: this
772) PetscInt :: chk_grp_id
773) type(option_type) :: option
774)
775) #if defined(SCORPIO_WRITE)
776) integer :: h5_chk_grp_id
777) integer, pointer :: dims(:)
778) integer, pointer :: start(:)
779) integer, pointer :: stride(:)
780) integer, pointer :: length(:)
781) integer :: timestepper_grp_id
782) #else
783) integer(HSIZE_T), pointer :: dims(:)
784) integer(HSIZE_T), pointer :: start(:)
785) integer(HSIZE_T), pointer :: stride(:)
786) integer(HSIZE_T), pointer :: length(:)
787) integer(HID_T) :: timestepper_grp_id
788) integer(HID_T) :: h5_chk_grp_id
789) #endif
790)
791) PetscMPIInt :: dataset_rank
792) character(len=MAXSTRINGLENGTH) :: dataset_name
793) character(len=MAXSTRINGLENGTH) :: string
794) PetscInt, pointer :: int_array(:)
795) PetscReal, pointer :: real_array(:)
796) PetscMPIInt :: hdf5_err
797)
798) string = "Timestepper"
799) h5_chk_grp_id = chk_grp_id
800) call h5gopen_f(h5_chk_grp_id, string, timestepper_grp_id, hdf5_err)
801)
802) allocate(start(1))
803) allocate(dims(1))
804) allocate(length(1))
805) allocate(stride(1))
806) allocate(int_array(1))
807) allocate(real_array(1))
808)
809) dataset_rank = 1
810) dims(1) = ONE_INTEGER
811) start(1) = 0
812) length(1) = ONE_INTEGER
813) stride(1) = ONE_INTEGER
814)
815) dataset_name = "Cumulative_newton_iterations" // CHAR(0)
816) call CheckPointReadIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
817) dims, start, length, stride, int_array, option)
818) this%cumulative_newton_iterations = int_array(1)
819)
820) dataset_name = "Cumulative_linear_iterations" // CHAR(0)
821) call CheckPointReadIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
822) dims, start, length, stride, int_array, option)
823) this%cumulative_linear_iterations = int_array(1)
824)
825) dataset_name = "Num_newton_iterations" // CHAR(0)
826) call CheckPointReadIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
827) dims, start, length, stride, int_array, option)
828) this%num_newton_iterations = int_array(1)
829)
830) dataset_name = "Time" // CHAR(0)
831) call CheckPointReadRealDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
832) dims, start, length, stride, real_array, option)
833) this%target_time = real_array(1)
834)
835) dataset_name = "Dt" // CHAR(0)
836) call CheckPointReadRealDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
837) dims, start, length, stride, real_array, option)
838) this%dt = real_array(1)
839)
840) dataset_name = "Prev_dt" // CHAR(0)
841) call CheckPointReadRealDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
842) dims, start, length, stride, real_array, option)
843) this%prev_dt = real_array(1)
844)
845) dataset_name = "Num_steps" // CHAR(0)
846) call CheckPointReadIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
847) dims, start, length, stride, int_array, option)
848) this%steps = int_array(1)
849)
850) dataset_name = "Cumulative_time_step_cuts" // CHAR(0)
851) call CheckPointReadIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
852) dims, start, length, stride, int_array, option)
853) this%cumulative_time_step_cuts = int_array(1)
854)
855) dataset_name = "Num_constant_time_steps" // CHAR(0)
856) call CheckPointReadIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
857) dims, start, length, stride, int_array, option)
858) this%num_constant_time_steps = int_array(1)
859)
860) dataset_name = "Num_contig_revert_due_to_sync" // CHAR(0)
861) call CheckPointReadIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
862) dims, start, length, stride, int_array, option)
863) this%num_contig_revert_due_to_sync = int_array(1)
864)
865) dataset_name = "Revert_dt" // CHAR(0)
866) call CheckPointReadIntDatasetHDF5(timestepper_grp_id, dataset_name, dataset_rank, &
867) dims, start, length, stride, int_array, option)
868) this%revert_dt = (int_array(1) == ONE_INTEGER)
869)
870) call h5gclose_f(timestepper_grp_id, hdf5_err)
871)
872) deallocate(start)
873) deallocate(dims)
874) deallocate(length)
875) deallocate(stride)
876) deallocate(int_array)
877) deallocate(real_array)
878)
879) end subroutine TimestepperBERestartHDF5
880) #endif
881)
882) ! ************************************************************************** !
883)
884) subroutine TimestepperBEGetHeader(this,header)
885) !
886) ! Gets values in checkpoint header.
887) !
888) ! Author: Glenn Hammond
889) ! Date: 07/25/13
890) !
891)
892) use Option_module
893)
894) implicit none
895)
896) #include "petsc/finclude/petscviewer.h"
897) #include "petsc/finclude/petscbag.h"
898)
899) class(timestepper_BE_type) :: this
900) class(stepper_BE_header_type) :: header
901)
902) this%cumulative_newton_iterations = header%cumulative_newton_iterations
903) this%cumulative_linear_iterations = header%cumulative_linear_iterations
904) this%num_newton_iterations = header%num_newton_iterations
905)
906) call TimestepperBaseGetHeader(this,header)
907)
908) end subroutine TimestepperBEGetHeader
909)
910) ! ************************************************************************** !
911)
912) subroutine TimestepperBEReset(this)
913) !
914) ! Zeros timestepper object members.
915) !
916) ! Author: Glenn Hammond
917) ! Date: 01/20/14
918) !
919)
920) implicit none
921)
922) class(timestepper_BE_type) :: this
923)
924) this%cumulative_newton_iterations = 0
925) this%cumulative_linear_iterations = 0
926) this%num_newton_iterations = 0
927)
928) call TimestepperBaseReset(this)
929)
930) end subroutine TimestepperBEReset
931)
932) ! ************************************************************************** !
933)
934) subroutine TimestepperBEPrintInfo(this,option)
935) !
936) ! Prints settings for base timestepper.
937) !
938) ! Author: Glenn Hammond
939) ! Date: 12/04/14
940) !
941) use Option_module
942)
943) implicit none
944)
945) class(timestepper_BE_type) :: this
946) type(option_type) :: option
947)
948) call TimestepperBasePrintInfo(this,option)
949) call SolverPrintNewtonInfo(this%solver,this%name,option)
950) call SolverPrintLinearInfo(this%solver,this%name,option)
951)
952) end subroutine TimestepperBEPrintInfo
953)
954) ! ************************************************************************** !
955)
956) subroutine TimestepperBEInputRecord(this)
957) !
958) ! Prints information about the time stepper to the input record.
959) ! To get a## format, must match that in simulation types.
960) !
961) ! Author: Jenn Frederick, SNL
962) ! Date: 03/17/2016
963) !
964)
965) implicit none
966)
967) class(timestepper_BE_type) :: this
968)
969) PetscInt :: id
970) character(len=MAXWORDLENGTH) :: word
971)
972) id = INPUT_RECORD_UNIT
973)
974) write(id,'(a29)',advance='no') 'pmc timestepper: '
975) write(id,'(a)') this%name
976)
977) write(id,'(a29)',advance='no') 'initial timestep size: '
978) write(word,*) this%dt_init
979) write(id,'(a)') trim(adjustl(word)) // ' sec'
980)
981) end subroutine TimestepperBEInputRecord
982)
983) ! ************************************************************************** !
984)
985) recursive subroutine TimestepperBEFinalizeRun(this,option)
986) !
987) ! Finalizes the time stepping
988) !
989) ! Author: Glenn Hammond
990) ! Date: 07/22/13
991) !
992)
993) use Option_module
994)
995) implicit none
996)
997) class(timestepper_BE_type) :: this
998) type(option_type) :: option
999)
1000) character(len=MAXSTRINGLENGTH) :: string
1001)
1002) #ifdef DEBUG
1003) call printMsg(option,'TimestepperBEFinalizeRun()')
1004) #endif
1005)
1006) if (OptionPrintToScreen(option)) then
1007) write(*,'(/,a," TS BE steps = ",i6," newton = ",i8," linear = ",i10, &
1008) & " cuts = ",i6)') &
1009) trim(this%name), &
1010) this%steps, &
1011) this%cumulative_newton_iterations, &
1012) this%cumulative_linear_iterations, &
1013) this%cumulative_time_step_cuts
1014) write(string,'(i12)') this%cumulative_wasted_linear_iterations
1015) write(*,'(a)') trim(this%name) // ' TS BE Wasted Linear Iterations = ' // &
1016) trim(adjustl(string))
1017) write(string,'(f12.1)') this%cumulative_solver_time
1018) write(*,'(a)') trim(this%name) // ' TS BE SNES time = ' // &
1019) trim(adjustl(string)) // ' seconds'
1020) endif
1021)
1022) end subroutine TimestepperBEFinalizeRun
1023)
1024) ! ************************************************************************** !
1025)
1026) subroutine TimestepperBEStrip(this)
1027) !
1028) ! Deallocates members of a time stepper
1029) !
1030) ! Author: Glenn Hammond
1031) ! Date: 07/22/13
1032) !
1033)
1034) implicit none
1035)
1036) class(timestepper_BE_type) :: this
1037)
1038) call TimestepperBaseStrip(this)
1039) call SolverDestroy(this%solver)
1040) call ConvergenceContextDestroy(this%convergence_context)
1041)
1042) if (associated(this%tfac)) deallocate(this%tfac)
1043) nullify(this%tfac)
1044)
1045) end subroutine TimestepperBEStrip
1046)
1047) ! ************************************************************************** !
1048)
1049) subroutine TimestepperBEDestroy(this)
1050) !
1051) ! Deallocates a time stepper
1052) !
1053) ! Author: Glenn Hammond
1054) ! Date: 07/22/13
1055) !
1056)
1057) implicit none
1058)
1059) class(timestepper_BE_type) :: this
1060)
1061) call TimestepperBEStrip(this)
1062)
1063) ! deallocate(this)
1064) ! nullify(this)
1065)
1066) end subroutine TimestepperBEDestroy
1067)
1068) end module Timestepper_BE_class