Octopus
states_mxll_restart.F90
Go to the documentation of this file.
1!! Copyright (C) 2019 R. Jestaedt, F. Bonafe, H. Appel, A. Rubio
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use batch_oct_m
25 use debug_oct_m
27 use global_oct_m
28 use grid_oct_m
31 use loct_oct_m
32 use mesh_oct_m
35 use mpi_oct_m
39 use parser_oct_m
43 use space_oct_m
47 use string_oct_m
48 use types_oct_m
50 use unit_oct_m
51
52 implicit none
53
54
55 private
56 public :: &
60
61contains
62
63 !----------------------------------------------------------
64 subroutine states_mxll_read_user_def(namespace, space, mesh, der, st, bc, user_def_rs_state)
65 type(namespace_t), intent(in) :: namespace
66 class(space_t), intent(in) :: space
67 class(mesh_t), intent(inout) :: mesh
68 type(derivatives_t), intent(in) :: der
69 type(states_mxll_t), intent(inout) :: st
70 type(bc_mxll_t), intent(inout) :: bc
71 complex(real64), contiguous, intent(inout) :: user_def_rs_state(:,:)
72
73 type(block_t) :: blk
74 integer :: il, nlines, idim, ncols, ip, state_from, ierr, maxwell_field
75 real(real64) :: xx(space%dim), rr, e_value, dummy, b_value
76 real(real64), allocatable :: e_field(:), b_field(:)
77 real(real64), allocatable :: total_efield(:,:), total_bfield(:,:)
78 complex(real64), allocatable :: rs_state_add(:), rs_state(:,:)
79 character(len=150), pointer :: filename_e_field, filename_b_field
80 character(1) :: cdim
81
82 integer, parameter :: &
83 STATE_FROM_FORMULA = 1, &
84 state_from_incident_waves = 2, &
85 state_from_file = -10010
86
88
89 call profiling_in('STATES_MXLL_READ_USER_DEF')
90
91 !%Variable UserDefinedInitialMaxwellStates
92 !%Type block
93 !%Section Maxwell
94 !%Description
95 !% The initial electromagnetic fields can be set by the user
96 !% with the <tt>UserDefinedMaxwellStates</tt> block variable.
97 !% The electromagnetic fields have to fulfill the
98 !% Maxwells equations in vacuum.
99 !%
100 !% Example:
101 !%
102 !% <tt>%UserDefinedMaxwellStates
103 !% <br>&nbsp;&nbsp; formula | 2 | "magnetic_field" | "-1/P_c * sin(x)"
104 !% <br>&nbsp;&nbsp; formula | 3 | "electric_field" | " sin(x) "
105 !% <br>%</tt>
106 !%
107 !% The second column specifies the component of the dimension of
108 !% the electric field and magnetic field. The first column
109 !% indicates that column four should be interpreted
110 !% as a formula for the corresponding state. P_c is the
111 !% speed of light constant.
112 !%
113 !% Alternatively, if column one states <tt>file</tt> the
114 !% electric field and magnetic field will be read from
115 !% the files given in column four.
116 !%
117 !% <tt>%UserDefinedMaxwellStates
118 !% <br>&nbsp;&nbsp; file | 3 | electric_field | "/path/to/file_electric_field_of_dimension_3"
119 !% <br>&nbsp;&nbsp; file | 2 | magnetic_field | "/path/to/file_magnetic_field_of_dimension_2"
120 !% <br>%</tt>
121 !%
122 !% The third option to define the initial state inside the box is to extend
123 !% the plane waves used as incident waves in the <tt>MaxwellIncidentWaves</tt> block,
124 !% as follows:
125 !%
126 !% <tt>%UserDefinedMaxwellStates
127 !% <br>&nbsp;&nbsp; use_incident_waves
128 !% <br>%</tt>
129 !%
130 !%Option file -10010
131 !% Read initial orbital from file.
132 !% Accepted file formats: obf, ncdf and csv.
133 !%Option formula 1
134 !% Calculate initial orbital by given analytic expression.
135 !%Option use_incident_waves 2
136 !% Extend the plane waves given in the <tt>MaxwellIncidentWaves</tt> block inside the box.
137 !%Option electric_field 1
138 !% This row defines the electric field component of the corresponding dimension
139 !%Option magnetic_field 2
140 !% This row defines the magnetic field component of the corresponding dimension
141 !%End
142
143 if (parse_block(namespace, 'UserDefinedInitialMaxwellStates', blk) == 0) then
144
145 safe_allocate(rs_state_add(1:mesh%np_part))
146 safe_allocate(rs_state(1:mesh%np, 1:3))
147
148 ! Set electromagnetic field equal to zero in the whole simulation box.
149 user_def_rs_state(:,:) = m_zero
150
151 ! find out how many lines (i.e. states) the block has
152 nlines = parse_block_n(blk)
153
154 write(message(1), '(a,i5)') 'Maxwell electromagnetic fields are added.'
155 write(message(2), '(a,i5)') ''
156 call messages_info(2, namespace=namespace)
157
158
159 safe_allocate(total_efield(1:mesh%np, 1:3))
160 safe_allocate(total_bfield(1:mesh%np, 1:3))
161 total_efield = m_zero
162 total_bfield = m_zero
163
164 ! read all lines
165 do il = 1, nlines
166 ! Check that number of columns is one or four.
167 ncols = parse_block_cols(blk, il - 1)
168 if (ncols /= 1 .and. ncols /= 4) then
169 message(1) = 'Each line in the UserDefinedMaxwellStates block must have'
170 message(2) = 'one or four columns.'
171 call messages_fatal(2, namespace=namespace)
172 end if
173
174 ! Calculate from expression, use or read from file?
175 call parse_block_integer(blk, il - 1, 0, state_from)
176
177 if (ncols > 1) then
178 call parse_block_integer(blk, il - 1, 1, idim)
179 write(cdim,'(I1)')idim
180 end if
181
182 rs_state_add(:) = m_zero
183
184 select case (state_from)
185
186 case (state_from_formula)
187
188 assert(ncols == 4)
189
190 ! parse formula string
191 call parse_block_integer(blk, il - 1, 2, maxwell_field)
192 if (maxwell_field == option__userdefinedinitialmaxwellstates__electric_field) then
193 call parse_block_string( blk, il - 1, 3, st%user_def_e_field(idim))
194 call messages_write(" E-field in dimension "//trim(cdim)//" : "//trim(st%user_def_e_field(idim)), fmt='(a,i1,2a)')
195 call conv_to_c_string(st%user_def_e_field(idim))
196 else if (maxwell_field == option__userdefinedinitialmaxwellstates__magnetic_field) then
197 call parse_block_string( blk, il - 1, 3, st%user_def_b_field(idim))
198 call messages_write(" B-field in dimension "//trim(cdim)//" : "//trim(st%user_def_b_field(idim)), fmt='(a,i1,2a)')
199 call conv_to_c_string(st%user_def_b_field(idim))
200 end if
201
202 if (maxwell_field == option__userdefinedinitialmaxwellstates__electric_field) then
203
204 ! fill Maxwell states with user-defined formulas
205 do ip = 1, mesh%np
206 xx = mesh%x(ip, :)
207 rr = norm2(xx)
208 ! parse user-defined expressions
209 call parse_expression(e_value, dummy, st%dim, xx, rr, m_zero, &
210 st%user_def_e_field(idim))
211 total_efield(ip, idim) = total_efield(ip, idim) + e_value
212 end do
213
214 else if (maxwell_field == option__userdefinedinitialmaxwellstates__magnetic_field) then
215 ! fill Maxwell states with user-defined formulas
216 do ip = 1, mesh%np
217 xx = mesh%x(ip, :)
218 rr = norm2(xx)
219
220 call parse_expression(b_value, dummy, st%dim, xx, rr, m_zero, &
221 st%user_def_b_field(idim))
222 total_bfield(ip, idim) = total_bfield(ip, idim) + b_value
223 end do
224
225 end if
226
227 case (state_from_file)
228
229 assert(ncols == 4)
230 ! The input format can be coded in column four now. As it is
231 ! not used now, we just say "file".
232 ! Read the filename.
233 call parse_block_integer(blk, il - 1, 2, maxwell_field)
234 safe_allocate(e_field(1:mesh%np))
235 safe_allocate(b_field(1:mesh%np))
236 e_field = m_zero
237 b_field = m_zero
238 if (maxwell_field == option__userdefinedinitialmaxwellstates__electric_field) then
239 call parse_block_string(blk, il - 1, 3, filename_e_field)
240 call messages_write(" E-field in dimension "//trim(cdim)//" : "//trim(filename_e_field), fmt='(a,i1,2a)')
241 call dio_function_input(filename_e_field, namespace, space, mesh, e_field(:), ierr)
242 if (ierr > 0) then
243 message(1) = 'Could not read the file!'
244 write(message(2),'(a,i1)') 'Error code: ', ierr
245 call messages_fatal(2, namespace=namespace)
246 end if
247 else if (maxwell_field == option__userdefinedinitialmaxwellstates__magnetic_field) then
248 call parse_block_string(blk, il - 1, 3, filename_b_field)
249 call messages_write(" B-field in dimension "//trim(cdim)//" : "//trim(filename_b_field), fmt='(a,i1,2a)')
250 call dio_function_input(filename_b_field, namespace, space, mesh, b_field(:), ierr)
251 if (ierr > 0) then
252 message(1) = 'Could not read the file!'
253 write(message(2),'(a,i1)') 'Error code: ', ierr
254 call messages_fatal(2, namespace=namespace)
255 end if
256 end if
257 ! fill state
258 call build_rs_vector(e_field(:), b_field(:), st%rs_sign, rs_state_add(:), &
259 st%ep(ip), st%mu(ip))
260
261 safe_deallocate_a(e_field)
262 safe_deallocate_a(b_field)
263
264 call lalg_axpy(mesh%np, m_one, rs_state_add, user_def_rs_state(:,idim))
265
266 case (state_from_incident_waves)
267
268 assert(ncols == 1)
269 call plane_waves_in_box_calculation(bc, m_zero, mesh, der, st, user_def_rs_state)
270
271 case default
272 message(1) = 'Wrong entry in UserDefinedMaxwellStates, column 2.'
273 message(2) = 'You may state "formula", "file" or "use_incident_waves" here.'
274 call messages_fatal(2, namespace=namespace)
275 end select
276
277 end do
278
279 if (state_from == state_from_formula) then
280 ! fill state
281 call build_rs_state(total_efield, total_bfield, st%rs_sign, rs_state, mesh, st%ep, st%mu)
282 call lalg_axpy(mesh%np, 3, m_one, rs_state, user_def_rs_state)
283 end if
284
285
286 safe_deallocate_a(total_efield)
287 safe_deallocate_a(total_bfield)
288
289 safe_deallocate_a(rs_state)
290 safe_deallocate_a(rs_state_add)
291 call parse_block_end(blk)
292 !call messages_print_with_emphasis(namespace=namespace)
293
294 else
295 call messages_variable_is_block(namespace, 'UserDefineInitialdStates')
296 end if
297
298 call profiling_out('STATES_MXLL_READ_USER_DEF')
299
301 end subroutine states_mxll_read_user_def
302
303
304 !----------------------------------------------------------
305 subroutine states_mxll_dump(restart, st, space, mesh, zff, zff_dim, ierr, iter, verbose)
306 type(restart_t), intent(in) :: restart
307 type(states_mxll_t), intent(in) :: st
308 class(space_t), intent(in) :: space
309 class(mesh_t), intent(in) :: mesh
310 complex(real64), intent(in) :: zff(:,:)
311 integer, intent(in) :: zff_dim
312 integer, intent(out) :: ierr
313 integer, optional, intent(in) :: iter
314 logical, optional, intent(in) :: verbose
315
316 integer :: iunit_wfns, iunit_states
317 integer :: err, err2(2), ist, idim, itot
318 integer :: root(1:p_strategy_max)
319 character(len=MAX_PATH_LEN) :: filename
320 character(len=300) :: lines(3)
321 logical :: should_write, verbose_
322
323 push_sub(states_mxll_dump)
324
325 verbose_ = optional_default(verbose, .true.)
326
327 ierr = 0
328
329 if (restart%skip()) then
330 pop_sub(states_mxll_dump)
331 return
332 end if
333
334 if (verbose_) then
335 message(1) = "Info: Writing Maxwell states."
336 call print_date(trim(message(1))//' ')
337 end if
338
339 call profiling_in("MAXWELL_RESTART_WRITE")
340
342
343 iunit_states = restart%open('maxwell_states')
344 write(lines(1),*) zff_dim
345 call restart%write(iunit_states, lines, 1, err)
346 if (err /= 0) ierr = ierr + 1
347 call restart%close(iunit_states)
348
349 iunit_wfns = restart%open('wfns')
350 lines(1) = '# #dim filename'
351 lines(2) = '%RS States'
352 call restart%write(iunit_wfns, lines, 2, err)
353 if (err /= 0) ierr = ierr + 2
354
355
356 itot = 1
357 root = 0
358 err2 = 0
359 ist = 1
360
361 do idim = 1, zff_dim
362 itot = itot + 1
363
364 root(p_strategy_domains) = mod(itot - 1, mesh%mpi_grp%size)
365 write(filename,'(i10.10)') itot
366
367 write(lines(1), '(i8,3a)') idim, ' | "', trim(filename), '"'
368 call restart%write(iunit_wfns, lines, 1, err)
369 if (err /= 0) err2(1) = err2(1) + 1
370
371 should_write = st%st_start <= ist .and. ist <= st%st_end
372
373 if (should_write) then
374 call restart%write_mesh_function(filename, mesh, zff(:,idim), err, root = root)
375 if (err /= 0) err2(2) = err2(2) + 1
376 end if
377 end do ! zff_dim
378
379 if (err2(1) /= 0) ierr = ierr + 8
380 if (err2(2) /= 0) ierr = ierr + 16
381
382 lines(1) = '%'
383 call restart%write(iunit_wfns, lines, 1, err)
384 if (err /= 0) ierr = ierr + 64
385 if (present(iter)) then
386 write(lines(1),'(a,i7)') 'Iter = ', iter
387 call restart%write(iunit_wfns, lines, 1, err)
388 if (err /= 0) ierr = ierr + 128
389 end if
390
391 call restart%close(iunit_wfns)
392
393 if (verbose_) then
394 message(1) = "Info: Finished writing Maxwell states."
395 call print_date(trim(message(1))//' ')
396 end if
397
399
400 call profiling_out("MAXWELL_RESTART_WRITE")
401 pop_sub(states_mxll_dump)
402 return
403
404 end subroutine states_mxll_dump
405
406 !----------------------------------------------------------
407 subroutine states_mxll_load(restart, st, mesh, namespace, space, zff, zff_dim, ierr, iter, lowest_missing, label, verbose)
408 type(restart_t), intent(in) :: restart
409 type(states_mxll_t), intent(inout) :: st
410 class(mesh_t), intent(in) :: mesh
411 type(namespace_t), intent(in) :: namespace
412 class(space_t), intent(in) :: space
413 complex(real64), contiguous, intent(inout) :: zff(:,:)
414 integer, intent(in) :: zff_dim
415 integer, intent(out) :: ierr
416 integer, optional, intent(out) :: iter
417 integer, optional, intent(out) :: lowest_missing(:)
418 character(len=*), optional, intent(in) :: label
419 logical, optional, intent(in) :: verbose
420
421 integer :: states_file, wfns_file, err, ist, idim, dim, mx_st_start, mx_st_end
422 integer :: idone, iread, ntodo
423 character(len=12) :: filename
424 character(len=1) :: char
425 logical, allocatable :: filled(:, :)
426 character(len=256) :: lines(3), label_
427
428 logical :: verbose_
429 character(len=256), allocatable :: restart_file(:, :)
430 logical, allocatable :: restart_file_present(:, :)
431
432 push_sub(states_mxll_load)
433
434 ierr = 0
435 dim = zff_dim
436 ist = 1
437
438 ! make sure these intent(out)`s are initialized no matter what
439 if (present(lowest_missing)) lowest_missing = 1
440 if (present(iter)) iter = 0
441
442 if (restart%skip()) then
443 ierr = -1
444 pop_sub(states_mxll_load)
445 return
446 end if
447
448 call profiling_in("RESTART_READ")
449
450 verbose_ = optional_default(verbose, .true.)
451
452 if (present(label)) then
453 label_ = trim(label)
454 end if
455
456 message(1) = 'Info: Reading Maxwell states'
457 if (len(trim(label_)) > 0) then
458 message(1) = trim(message(1)) // trim(label_)
459 end if
460 message(1) = trim(message(1)) // "."
461 if (verbose_) call print_date(trim(message(1))//' ')
462
463 states_file = restart%open('maxwell_states')
464 call restart%read(states_file, lines, 1, err)
465 if (err /= 0) then
466 ierr = ierr - 2
467 else
468 read(lines(1), *) idim
469 end if
470 call restart%close(states_file)
471
472 ! open files to read
473 wfns_file = restart%open('wfns')
474 call restart%read(wfns_file, lines, 2, err)
475 if (err /= 0) then
476 ierr = ierr - 2**5
477 end if
478
479 ! If any error occured up to this point then it is not worth continuing,
480 ! as there something fundamentally wrong with the restart files
481 if (ierr /= 0) then
482 call restart%close(wfns_file)
483 call profiling_out("RESTART_READ")
484 pop_sub(states_mxll_load)
485 return
486 end if
487
488 safe_allocate(restart_file(1:zff_dim, st%st_start:st%st_end))
489 safe_allocate(restart_file_present(1:zff_dim,st%st_start:st%st_end))
490 restart_file_present = .false.
491
492 ! Next we read the list of states from the files.
493 ! Errors in reading the information of a specific state from the files are ignored
494 ! at this point, because later we will skip reading the wavefunction of that state.
495 do
496 call restart%read(wfns_file, lines, 1, err)
497 if (err == 0) then
498 read(lines(1), '(a)') char
499 if (char == '%') then
500 !We reached the end of the file
501 exit
502 else
503 read(lines(1), *) idim, char, filename
504 end if
505 end if
506
507 if (ist >= st%st_start .and. ist <= st%st_end) then
508 restart_file(idim, ist) = trim(filename)
509 restart_file_present(idim, ist) = .true.
510 end if
511
512 end do
513 if (present(iter)) then
514 call restart%read(wfns_file, lines, 1, err)
515 if (err /= 0) then
516 ierr = ierr - 2**8
517 else
518 read(lines(1), *) filename, filename, iter
519 end if
520 end if
521
522 call restart%close(wfns_file)
523
524 ! Now we read the wavefunctions. At this point we need to have all the information from the
525 ! states, and wfns files in order to avoid serialisation of reads, as restart_read
526 ! works as a barrier.
527
528 mx_st_start=st%st_start
529 mx_st_end=st%st_end
530 safe_allocate(filled(1:zff_dim,mx_st_start:mx_st_end))
531 filled = .false.
532
533 if (present(lowest_missing)) lowest_missing = st%nst + 1
534
535 iread = 0
536 if (mpi_world%is_root() .and. verbose_) then
537 idone = 0
538 ntodo = st%lnst*zff_dim
539 call loct_progress_bar(-1, ntodo)
540 end if
541
542 ist = 1
543 do idim = 1, zff_dim
544 if (.not. restart_file_present(idim, ist)) then
545 if (present(lowest_missing)) then
546 lowest_missing(idim) = min(lowest_missing(idim), ist)
547 end if
548 cycle
549 end if
550
551 call restart%read_mesh_function(restart_file(idim, ist), mesh, zff(:,idim), err)
552
553 if (err == 0) then
554 filled(idim, ist) = .true.
555 iread = iread + 1
556 else if (present(lowest_missing)) then
557 lowest_missing(idim) = min(lowest_missing(idim), ist)
558 end if
559
560 if (mpi_world%is_root() .and. verbose_) then
561 idone = idone + 1
562 call loct_progress_bar(idone, ntodo)
563 end if
564
565 end do
566
567 safe_deallocate_a(restart_file)
568 safe_deallocate_a(restart_file_present)
569 safe_deallocate_a(filled)
570
571 if (mpi_world%is_root() .and. verbose_) then
572 call messages_new_line()
573 end if
574
575 if (ierr == 0 .and. iread /= st%nst * zff_dim) then
576 if (iread > 0) then
577 ierr = iread
578 else
579 ierr = -1
580 end if
581 ! otherwise ierr = 0 would mean either all was read correctly, or nothing at all was read!
582
583 call messages_print_with_emphasis(msg='Reading Maxwell states.', namespace=namespace)
584 write(message(1),'(a,i6,a,i6,a)') 'Only ', iread,' files out of ', &
585 st%nst * zff_dim, ' could be read.'
586 call messages_info(1, namespace=namespace)
587 call messages_print_with_emphasis(namespace=namespace)
588 end if
589
590 message(1) = 'Info: Maxwell states reading done.'
591 if (verbose_) call print_date(trim(message(1))//' ')
592
593 call profiling_out("RESTART_READ")
594 pop_sub(states_mxll_load)
595
596 end subroutine states_mxll_load
597
599
600
601!! Local Variables:
602!! mode: f90
603!! coding: utf-8
604!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:173
block signals while writing the restart files
Definition: restart.F90:318
unblock signals when writing restart is finished
Definition: restart.F90:325
This module implements batches of mesh functions.
Definition: batch.F90:135
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
This module handles the calculation mode.
integer, parameter, public p_strategy_max
integer, parameter, public p_strategy_domains
parallelization in domains
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_one
Definition: global.F90:192
This module implements the underlying real-space grid.
Definition: grid.F90:119
subroutine, public dio_function_input(filename, namespace, space, mesh, ff, ierr, map)
Reads a mesh function from file filename, and puts it into ff. If the map argument is passed,...
System information (time, memory, sysname)
Definition: loct.F90:117
subroutine, public loct_progress_bar(a, maxcount)
A wrapper around the progress bar, such that it can be silenced without needing to dress the call wit...
Definition: loct.F90:276
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
character(len=512), private msg
Definition: messages.F90:167
subroutine, public messages_variable_is_block(namespace, name)
Definition: messages.F90:1047
subroutine, public print_date(str)
Definition: messages.F90:983
subroutine, public messages_new_line()
Definition: messages.F90:1112
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
subroutine, public parse_block_string(blk, l, c, res, convert_to_c)
Definition: parser.F90:810
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:615
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
subroutine, public plane_waves_in_box_calculation(bc, time, mesh, der, st, rs_state)
This module handles spin dimensions of the states and the k-point distribution.
This module handles reading and writing restart information for the states_elec_t.
subroutine, public build_rs_vector(e_vector, b_vector, rs_sign, rs_vector, ep_element, mu_element)
subroutine, public build_rs_state(e_field, b_field, rs_sign, rs_state, mesh, ep_field, mu_field, np)
subroutine, public states_mxll_load(restart, st, mesh, namespace, space, zff, zff_dim, ierr, iter, lowest_missing, label, verbose)
subroutine, public states_mxll_dump(restart, st, space, mesh, zff, zff_dim, ierr, iter, verbose)
subroutine, public states_mxll_read_user_def(namespace, space, mesh, der, st, bc, user_def_rs_state)
subroutine, public conv_to_c_string(str)
converts to c string
Definition: string.F90:257
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
Describes mesh distribution to nodes.
Definition: mesh.F90:187
int true(void)