rungvec.f90 Source File


This file depends on

sourcefile~~rungvec.f90~~EfferentGraph sourcefile~rungvec.f90 rungvec.f90 sourcefile~analyze.f90 analyze.f90 sourcefile~rungvec.f90->sourcefile~analyze.f90 sourcefile~functional.f90 functional.f90 sourcefile~rungvec.f90->sourcefile~functional.f90 sourcefile~globals.f90 globals.f90 sourcefile~rungvec.f90->sourcefile~globals.f90 sourcefile~mhd3d_vars.f90 mhd3d_vars.f90 sourcefile~rungvec.f90->sourcefile~mhd3d_vars.f90 sourcefile~output.f90 output.f90 sourcefile~rungvec.f90->sourcefile~output.f90 sourcefile~readintools.f90 readintools.f90 sourcefile~rungvec.f90->sourcefile~readintools.f90 sourcefile~restart.f90 restart.f90 sourcefile~rungvec.f90->sourcefile~restart.f90 sourcefile~analyze.f90->sourcefile~globals.f90 sourcefile~analyze.f90->sourcefile~mhd3d_vars.f90 sourcefile~analyze.f90->sourcefile~readintools.f90 sourcefile~analyze_vars.f90 analyze_vars.f90 sourcefile~analyze.f90->sourcefile~analyze_vars.f90 sourcefile~cubic_spline.f90 cubic_spline.f90 sourcefile~analyze.f90->sourcefile~cubic_spline.f90 sourcefile~mod_mpi.f90 mod_mpi.f90 sourcefile~analyze.f90->sourcefile~mod_mpi.f90 sourcefile~output_csv.f90 output_csv.f90 sourcefile~analyze.f90->sourcefile~output_csv.f90 sourcefile~output_vars.f90 output_vars.f90 sourcefile~analyze.f90->sourcefile~output_vars.f90 sourcefile~output_vtk.f90 output_vtk.f90 sourcefile~analyze.f90->sourcefile~output_vtk.f90 sourcefile~vmec.f90 vmec.f90 sourcefile~analyze.f90->sourcefile~vmec.f90 sourcefile~vmec_vars.f90 vmec_vars.f90 sourcefile~analyze.f90->sourcefile~vmec_vars.f90 sourcefile~write_modes.f90 write_modes.f90 sourcefile~analyze.f90->sourcefile~write_modes.f90 sourcefile~functional.f90->sourcefile~globals.f90 sourcefile~c_functional.f90 c_functional.f90 sourcefile~functional.f90->sourcefile~c_functional.f90 sourcefile~mhd3d.f90 mhd3d.f90 sourcefile~functional.f90->sourcefile~mhd3d.f90 sourcefile~mhd3d_vars.f90->sourcefile~globals.f90 sourcefile~base.f90 base.f90 sourcefile~mhd3d_vars.f90->sourcefile~base.f90 sourcefile~c_rprofile.f90 c_rprofile.f90 sourcefile~mhd3d_vars.f90->sourcefile~c_rprofile.f90 sourcefile~hmap.f90 hmap.f90 sourcefile~mhd3d_vars.f90->sourcefile~hmap.f90 sourcefile~sgrid.f90 sgrid.f90 sourcefile~mhd3d_vars.f90->sourcefile~sgrid.f90 sourcefile~sol_var_mhd3d.f90 sol_var_mhd3d.f90 sourcefile~mhd3d_vars.f90->sourcefile~sol_var_mhd3d.f90 sourcefile~output.f90->sourcefile~globals.f90 sourcefile~output.f90->sourcefile~readintools.f90 sourcefile~output.f90->sourcefile~output_vars.f90 sourcefile~readintools.f90->sourcefile~globals.f90 sourcefile~readintools.f90->sourcefile~mod_mpi.f90 sourcefile~restart.f90->sourcefile~globals.f90 sourcefile~restart.f90->sourcefile~mhd3d_vars.f90 sourcefile~restart.f90->sourcefile~base.f90 sourcefile~mhd3d_evalfunc.f90 mhd3d_evalfunc.f90 sourcefile~restart.f90->sourcefile~mhd3d_evalfunc.f90 sourcefile~restart.f90->sourcefile~mod_mpi.f90 sourcefile~restart.f90->sourcefile~output_vars.f90 sourcefile~readstate.f90 readstate.f90 sourcefile~restart.f90->sourcefile~readstate.f90 sourcefile~readstate_vars.f90 readstate_vars.f90 sourcefile~restart.f90->sourcefile~readstate_vars.f90 sourcefile~restart_vars.f90 restart_vars.f90 sourcefile~restart.f90->sourcefile~restart_vars.f90 sourcefile~restart.f90->sourcefile~sgrid.f90 sourcefile~restart.f90->sourcefile~sol_var_mhd3d.f90 sourcefile~analyze_vars.f90->sourcefile~globals.f90 sourcefile~base.f90->sourcefile~globals.f90 sourcefile~base.f90->sourcefile~sgrid.f90 sourcefile~fbase.f90 fbase.f90 sourcefile~base.f90->sourcefile~fbase.f90 sourcefile~sbase.f90 sbase.f90 sourcefile~base.f90->sourcefile~sbase.f90 sourcefile~c_functional.f90->sourcefile~globals.f90 sourcefile~c_rprofile.f90->sourcefile~globals.f90 sourcefile~cubic_spline.f90->sourcefile~globals.f90 sourcefile~sll_m_bsplines.f90 sll_m_bsplines.f90 sourcefile~cubic_spline.f90->sourcefile~sll_m_bsplines.f90 sourcefile~sll_m_spline_matrix.f90 sll_m_spline_matrix.f90 sourcefile~cubic_spline.f90->sourcefile~sll_m_spline_matrix.f90 sourcefile~hmap.f90->sourcefile~globals.f90 sourcefile~c_hmap.f90 c_hmap.f90 sourcefile~hmap.f90->sourcefile~c_hmap.f90 sourcefile~hmap_axisnb.f90 hmap_axisNB.f90 sourcefile~hmap.f90->sourcefile~hmap_axisnb.f90 sourcefile~hmap_cyl.f90 hmap_cyl.f90 sourcefile~hmap.f90->sourcefile~hmap_cyl.f90 sourcefile~hmap_frenet.f90 hmap_frenet.f90 sourcefile~hmap.f90->sourcefile~hmap_frenet.f90 sourcefile~hmap_knot.f90 hmap_knot.f90 sourcefile~hmap.f90->sourcefile~hmap_knot.f90 sourcefile~hmap_rz.f90 hmap_RZ.f90 sourcefile~hmap.f90->sourcefile~hmap_rz.f90 sourcefile~mhd3d.f90->sourcefile~analyze.f90 sourcefile~mhd3d.f90->sourcefile~globals.f90 sourcefile~mhd3d.f90->sourcefile~mhd3d_vars.f90 sourcefile~mhd3d.f90->sourcefile~readintools.f90 sourcefile~mhd3d.f90->sourcefile~restart.f90 sourcefile~mhd3d.f90->sourcefile~base.f90 sourcefile~mhd3d.f90->sourcefile~c_functional.f90 sourcefile~mhd3d.f90->sourcefile~c_rprofile.f90 sourcefile~mhd3d.f90->sourcefile~cubic_spline.f90 sourcefile~mhd3d.f90->sourcefile~hmap.f90 sourcefile~mhd3d.f90->sourcefile~mhd3d_evalfunc.f90 sourcefile~mhd3d.f90->sourcefile~mod_mpi.f90 sourcefile~mhd3d.f90->sourcefile~output_vars.f90 sourcefile~mhd3d.f90->sourcefile~restart_vars.f90 sourcefile~mhd3d.f90->sourcefile~sgrid.f90 sourcefile~mhd3d.f90->sourcefile~sol_var_mhd3d.f90 sourcefile~mhd3d.f90->sourcefile~vmec.f90 sourcefile~mhd3d.f90->sourcefile~vmec_vars.f90 sourcefile~boundaryfromfile.f90 boundaryFromFile.f90 sourcefile~mhd3d.f90->sourcefile~boundaryfromfile.f90 sourcefile~lambda_solve.f90 lambda_solve.f90 sourcefile~mhd3d.f90->sourcefile~lambda_solve.f90 sourcefile~rprofile_bspline.f90 rprofile_bspline.f90 sourcefile~mhd3d.f90->sourcefile~rprofile_bspline.f90 sourcefile~rprofile_polynomial.f90 rprofile_polynomial.f90 sourcefile~mhd3d.f90->sourcefile~rprofile_polynomial.f90 sourcefile~mhd3d_evalfunc.f90->sourcefile~globals.f90 sourcefile~mhd3d_evalfunc.f90->sourcefile~mhd3d_vars.f90 sourcefile~mhd3d_evalfunc.f90->sourcefile~base.f90 sourcefile~mhd3d_evalfunc.f90->sourcefile~mod_mpi.f90 sourcefile~mhd3d_evalfunc.f90->sourcefile~sol_var_mhd3d.f90 sourcefile~mhd3d_evalfunc.f90->sourcefile~sll_m_spline_matrix.f90 sourcefile~sll_m_spline_matrix_banded.f90 sll_m_spline_matrix_banded.f90 sourcefile~mhd3d_evalfunc.f90->sourcefile~sll_m_spline_matrix_banded.f90 sourcefile~mod_mpi.f90->sourcefile~globals.f90 sourcefile~output_csv.f90->sourcefile~globals.f90 sourcefile~output_vars.f90->sourcefile~globals.f90 sourcefile~output_vtk.f90->sourcefile~globals.f90 sourcefile~readstate.f90->sourcefile~globals.f90 sourcefile~readstate.f90->sourcefile~base.f90 sourcefile~readstate.f90->sourcefile~hmap.f90 sourcefile~readstate.f90->sourcefile~readstate_vars.f90 sourcefile~readstate.f90->sourcefile~sgrid.f90 sourcefile~readstate.f90->sourcefile~fbase.f90 sourcefile~readstate.f90->sourcefile~sbase.f90 sourcefile~readstate_vars.f90->sourcefile~globals.f90 sourcefile~readstate_vars.f90->sourcefile~base.f90 sourcefile~readstate_vars.f90->sourcefile~hmap.f90 sourcefile~readstate_vars.f90->sourcefile~sgrid.f90 sourcefile~readstate_vars.f90->sourcefile~sbase.f90 sourcefile~restart_vars.f90->sourcefile~globals.f90 sourcefile~sgrid.f90->sourcefile~globals.f90 sourcefile~sol_var_mhd3d.f90->sourcefile~globals.f90 sourcefile~c_sol_var.f90 c_sol_var.f90 sourcefile~sol_var_mhd3d.f90->sourcefile~c_sol_var.f90 sourcefile~vmec.f90->sourcefile~globals.f90 sourcefile~vmec.f90->sourcefile~readintools.f90 sourcefile~vmec.f90->sourcefile~cubic_spline.f90 sourcefile~vmec.f90->sourcefile~vmec_vars.f90 sourcefile~vmec.f90->sourcefile~rprofile_bspline.f90 sourcefile~vmec_vars.f90->sourcefile~globals.f90 sourcefile~vmec_vars.f90->sourcefile~c_rprofile.f90 sourcefile~vmec_vars.f90->sourcefile~cubic_spline.f90 sourcefile~write_modes.f90->sourcefile~globals.f90 sourcefile~write_modes.f90->sourcefile~output_csv.f90 sourcefile~boundaryfromfile.f90->sourcefile~globals.f90 sourcefile~boundaryfromfile.f90->sourcefile~fbase.f90 sourcefile~io_netcdf.f90 io_netcdf.f90 sourcefile~boundaryfromfile.f90->sourcefile~io_netcdf.f90 sourcefile~c_hmap.f90->sourcefile~globals.f90 sourcefile~c_sol_var.f90->sourcefile~globals.f90 sourcefile~fbase.f90->sourcefile~globals.f90 sourcefile~hmap_axisnb.f90->sourcefile~globals.f90 sourcefile~hmap_axisnb.f90->sourcefile~readintools.f90 sourcefile~hmap_axisnb.f90->sourcefile~analyze_vars.f90 sourcefile~hmap_axisnb.f90->sourcefile~mod_mpi.f90 sourcefile~hmap_axisnb.f90->sourcefile~output_csv.f90 sourcefile~hmap_axisnb.f90->sourcefile~output_vtk.f90 sourcefile~hmap_axisnb.f90->sourcefile~c_hmap.f90 sourcefile~hmap_axisnb.f90->sourcefile~fbase.f90 sourcefile~hmap_axisnb.f90->sourcefile~io_netcdf.f90 sourcefile~output_netcdf.f90 output_netcdf.f90 sourcefile~hmap_axisnb.f90->sourcefile~output_netcdf.f90 sourcefile~hmap_cyl.f90->sourcefile~globals.f90 sourcefile~hmap_cyl.f90->sourcefile~readintools.f90 sourcefile~hmap_cyl.f90->sourcefile~c_hmap.f90 sourcefile~hmap_frenet.f90->sourcefile~globals.f90 sourcefile~hmap_frenet.f90->sourcefile~readintools.f90 sourcefile~hmap_frenet.f90->sourcefile~analyze_vars.f90 sourcefile~hmap_frenet.f90->sourcefile~output_csv.f90 sourcefile~hmap_frenet.f90->sourcefile~output_vtk.f90 sourcefile~hmap_frenet.f90->sourcefile~c_hmap.f90 sourcefile~hmap_frenet.f90->sourcefile~output_netcdf.f90 sourcefile~hmap_knot.f90->sourcefile~globals.f90 sourcefile~hmap_knot.f90->sourcefile~readintools.f90 sourcefile~hmap_knot.f90->sourcefile~c_hmap.f90 sourcefile~hmap_rz.f90->sourcefile~globals.f90 sourcefile~hmap_rz.f90->sourcefile~c_hmap.f90 sourcefile~lambda_solve.f90->sourcefile~globals.f90 sourcefile~lambda_solve.f90->sourcefile~base.f90 sourcefile~lambda_solve.f90->sourcefile~hmap.f90 sourcefile~lambda_solve.f90->sourcefile~fbase.f90 sourcefile~linalg.f90 linalg.f90 sourcefile~lambda_solve.f90->sourcefile~linalg.f90 sourcefile~rprofile_bspline.f90->sourcefile~globals.f90 sourcefile~rprofile_bspline.f90->sourcefile~c_rprofile.f90 sourcefile~rprofile_bspline.f90->sourcefile~sll_m_bsplines.f90 sourcefile~rprofile_polynomial.f90->sourcefile~globals.f90 sourcefile~rprofile_polynomial.f90->sourcefile~c_rprofile.f90 sourcefile~sbase.f90->sourcefile~globals.f90 sourcefile~sbase.f90->sourcefile~sgrid.f90 sourcefile~sbase.f90->sourcefile~sll_m_bsplines.f90 sourcefile~sbase.f90->sourcefile~sll_m_spline_matrix.f90 sourcefile~basis1d.f90 basis1d.f90 sourcefile~sbase.f90->sourcefile~basis1d.f90 sourcefile~sbase.f90->sourcefile~linalg.f90 sourcefile~sll_m_boundary_condition_descriptors.f90 sll_m_boundary_condition_descriptors.f90 sourcefile~sbase.f90->sourcefile~sll_m_boundary_condition_descriptors.f90 sourcefile~sll_m_spline_interpolator_1d.f90 sll_m_spline_interpolator_1d.f90 sourcefile~sbase.f90->sourcefile~sll_m_spline_interpolator_1d.f90 sourcefile~sll_m_bsplines_base.f90 sll_m_bsplines_base.f90 sourcefile~sll_m_bsplines.f90->sourcefile~sll_m_bsplines_base.f90 sourcefile~sll_m_bsplines_non_uniform.f90 sll_m_bsplines_non_uniform.f90 sourcefile~sll_m_bsplines.f90->sourcefile~sll_m_bsplines_non_uniform.f90 sourcefile~sll_m_bsplines_uniform.f90 sll_m_bsplines_uniform.f90 sourcefile~sll_m_bsplines.f90->sourcefile~sll_m_bsplines_uniform.f90 sourcefile~sll_m_working_precision.f90 sll_m_working_precision.f90 sourcefile~sll_m_bsplines.f90->sourcefile~sll_m_working_precision.f90 sourcefile~sll_m_spline_matrix.f90->sourcefile~sll_m_spline_matrix_banded.f90 sourcefile~sll_m_spline_matrix_base.f90 sll_m_spline_matrix_base.f90 sourcefile~sll_m_spline_matrix.f90->sourcefile~sll_m_spline_matrix_base.f90 sourcefile~sll_m_spline_matrix_dense.f90 sll_m_spline_matrix_dense.f90 sourcefile~sll_m_spline_matrix.f90->sourcefile~sll_m_spline_matrix_dense.f90 sourcefile~sll_m_spline_matrix.f90->sourcefile~sll_m_working_precision.f90 sourcefile~sll_m_spline_matrix_banded.f90->sourcefile~sll_m_spline_matrix_base.f90 sourcefile~sll_m_spline_matrix_banded.f90->sourcefile~sll_m_working_precision.f90 sourcefile~basis1d.f90->sourcefile~globals.f90 sourcefile~basis1d.f90->sourcefile~linalg.f90 sourcefile~io_netcdf.f90->sourcefile~globals.f90 sourcefile~linalg.f90->sourcefile~globals.f90 sourcefile~output_netcdf.f90->sourcefile~globals.f90 sourcefile~output_netcdf.f90->sourcefile~io_netcdf.f90 sourcefile~sll_m_bsplines_base.f90->sourcefile~sll_m_working_precision.f90 sourcefile~sll_m_bsplines_non_uniform.f90->sourcefile~sll_m_bsplines_base.f90 sourcefile~sll_m_bsplines_non_uniform.f90->sourcefile~sll_m_working_precision.f90 sourcefile~sll_m_bsplines_uniform.f90->sourcefile~sll_m_bsplines_base.f90 sourcefile~sll_m_bsplines_uniform.f90->sourcefile~sll_m_working_precision.f90 sourcefile~sll_m_spline_interpolator_1d.f90->sourcefile~sll_m_spline_matrix.f90 sourcefile~sll_m_spline_interpolator_1d.f90->sourcefile~sll_m_boundary_condition_descriptors.f90 sourcefile~sll_m_spline_interpolator_1d.f90->sourcefile~sll_m_bsplines_base.f90 sourcefile~sll_m_spline_interpolator_1d.f90->sourcefile~sll_m_working_precision.f90 sourcefile~sll_m_spline_1d.f90 sll_m_spline_1d.f90 sourcefile~sll_m_spline_interpolator_1d.f90->sourcefile~sll_m_spline_1d.f90 sourcefile~sll_m_spline_matrix_base.f90->sourcefile~sll_m_working_precision.f90 sourcefile~sll_m_spline_matrix_dense.f90->sourcefile~sll_m_spline_matrix_base.f90 sourcefile~sll_m_spline_matrix_dense.f90->sourcefile~sll_m_working_precision.f90 sourcefile~sll_m_spline_1d.f90->sourcefile~sll_m_bsplines_base.f90 sourcefile~sll_m_spline_1d.f90->sourcefile~sll_m_working_precision.f90

Files dependent on this one

sourcefile~~rungvec.f90~~AfferentGraph sourcefile~rungvec.f90 rungvec.f90 sourcefile~gvec.f90 gvec.f90 sourcefile~gvec.f90->sourcefile~rungvec.f90 sourcefile~run.f90 run.f90 sourcefile~run.f90->sourcefile~rungvec.f90

Source Code

!===================================================================================================================================
! Copyright (c) 2025 GVEC Contributors, Max Planck Institute for Plasma Physics
! License: MIT
!===================================================================================================================================
#include "defines.h"


!===================================================================================================================================
!>
!!# **GVEC** Driver Module
!!
!===================================================================================================================================
MODULE MODgvec_rungvec
IMPLICIT NONE
PUBLIC

!INTERFACE rungvec
!  MODULE PROCEDURE rungvec
!END INTERFACE
!===================================================================================================================================
CONTAINS

SUBROUTINE rungvec(parameterFile,restartfile_in)
USE_MPI
USE MODgvec_Globals    ,ONLY : wp,fmt_sep,n_warnings_occured,testdbg,testlevel,testUnit,nfailedMsg,UNIT_stdOut,GETFREEUNIT
USE MODgvec_Globals    ,ONLY : GetTime,MPIRoot,nRanks,myRank
USE MODgvec_Analyze    ,ONLY : InitAnalyze,FinalizeAnalyze
USE MODgvec_Output     ,ONLY : InitOutput,FinalizeOutput
USE MODgvec_Restart    ,ONLY : InitRestart,FinalizeRestart
USE MODgvec_ReadInTools,ONLY : FillStrings,GETLOGICAL,GETINT,IgnoredStrings,FinalizeReadIn
USE MODgvec_Functional ,ONLY : t_functional, InitFunctional,FinalizeFunctional
USE MODgvec_MHD3D_Vars ,ONLY : maxIter
!$ USE omp_lib
IMPLICIT NONE
!-----------------------------------------------------------------------------------------------------------------------------------
!INPUT VARIABLES
  CHARACTER(LEN=*),INTENT(IN)             :: Parameterfile  !! input parameters of GVEC
  CHARACTER(LEN=*),INTENT(IN),OPTIONAL    :: RestartFile_in  !! if present, a restart will be executed
!-----------------------------------------------------------------------------------------------------------------------------------
!local variables
INTEGER                 :: which_functional
INTEGER                 :: TimeArray(8)
CHARACTER(LEN=255)      :: testfile
REAL(wp)                :: StartTimeTotal,EndTimeTotal,StartTime,EndTime
CLASS(t_functional),ALLOCATABLE   :: functional
LOGICAL                 :: dorestart
!===================================================================================================================================
  __PERFINIT
  __PERFON('main')


  StartTimeTotal=GetTime()
  SWRITE(Unit_stdOut,fmt_sep)
  CALL DATE_AND_TIME(values=TimeArray) ! get System time
  SWRITE(UNIT_stdOut,'(A,I4.2,"-",I2.2,"-",I2.2,1X,I2.2,":",I2.2,":",I2.2)') &
    '%%% Sys date : ',timeArray(1:3),timeArray(5:7)

  !header
  SWRITE(Unit_stdOut,fmt_sep)
  SWRITE(Unit_stdOut,'(18(("*",A128,2X,"*",:,"\n")))')&
 '  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - '&
,' - - - - - - - - - - - - GGGGGGGGGGGGGGG - VVVVVVVV  - - - -  VVVVVVVV - EEEEEEEEEEEEEEEEEEEEEE  - - - - CCCCCCCCCCCCCCC - -  '&
,'  - - - - - - - - - - GGGG::::::::::::G - V::::::V  - - - -  V::::::V - E::::::::::::::::::::E  - - - CCCC::::::::::::C - - - '&
,' - - - - - - - - - GGG:::::::::::::::G - V::::::V  - - - -  V::::::V - E::::::::::::::::::::E  - - CCC:::::::::::::::C - - -  '&
,'  - - - - - - -  GG:::::GGGGGGGG::::G - V::::::V  - - - -  V::::::V - EEE:::::EEEEEEEEE::::E  -  CC:::::CCCCCCCC::::C - - - - '&
,' - - - - - - - GG:::::GG  - - GGGGGG -  V:::::V  - - - -  V:::::V  - - E:::::E - - - EEEEEE  - CC:::::CC - -  CCCCCC - - - -  '&
,'  - - - - - - G:::::GG  - - - - - - - - V:::::V - - - - V:::::V - - - E:::::E - - - - - - - - C:::::CC  - - - - - - - - - - - '&
,' - - - - - - G:::::G - - - - - - - - -  V:::::V  - -  V:::::V  - - - E:::::EEEEEEEEEEE - - - C:::::C - - - - - - - - - - - -  '&
,'  - - - - - G:::::G -  GGGGGGGGGG - - - V:::::V - - V:::::V - - - - E:::::::::::::::E - - - C:::::C - - - - - - - - - - - - - '&
,' - - - - - G:::::G -  G::::::::G - - -  V:::::V   V:::::V  - - - - E:::::::::::::::E - - - C:::::C - - - - - - - - - - - - -  '&
,'  - - - - G:::::G -  GGGGG::::G - - - - V:::::V V:::::V - - - - - E:::::EEEEEEEEEEE - - - C:::::C - - - - - - - - - - - - - - '&
,' - - - - G:::::G - - -  G::::G - - - -  V:::::V:::::V  - - - - - E:::::E - - - - - - - - C:::::C - - - - - - - - - - - - - -  '&
,'  - - -  G:::::G  - -  G::::G - - - - - V:::::::::V - - - - - - E:::::E - - - EEEEEE  -  C:::::C  - -  CCCCCC - - - - - - - - '&
,' - - - - G::::::GGGGGGG::::G - - - - -  V:::::::V  - - - - - EEE:::::EEEEEEEEE::::E  - - C::::::CCCCCCC::::C - - - - - - - -  '&
,'  - - - - G:::::::::::::::G - - - - - - V:::::V - - - - - - E::::::::::::::::::::E  - - - C:::::::::::::::C - - - - - - - - - '&
,' - - - - - GG::::GGGG::::G - - - - - -  V:::V  - - - - - - E::::::::::::::::::::E  - - - - CC::::::::::::C - - - - - - - - -  '&
,'  - - - - -  GGGG  GGGGGG - - - - - - - VVV - - - - - - - EEEEEEEEEEEEEEEEEEEEEE  - - - - -  CCCCCCCCCCCC - - - - - - - - - - '&
,' - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -  '
  SWRITE(Unit_stdOut,fmt_sep)
  !.only executes if compiled with OpenMP
!$ SWRITE(UNIT_stdOut,'(A,I6)')'   Number of OpenMP threads : ',OMP_GET_MAX_THREADS()
!$ SWRITE(Unit_stdOut,'(132("="))')
  !.only executes if compiled with MPI
# if MPI
  SWRITE(UNIT_stdOut,'(A,I6)')'   Number of MPI tasks : ',nRanks
  SWRITE(Unit_stdOut,fmt_sep)
# endif
#include  "configuration-cmake.f90"
  SWRITE(Unit_stdOut,fmt_sep)

  CALL FillStrings(ParameterFile) !< readin parameterfile, done on MPI root + Bcast

  testdbg =GETLOGICAL('testdbg',Proposal=.FALSE.)
  testlevel=GETINT('testlevel',Proposal=-1)
  IF(testlevel.GT.0)THEN
    testUnit=GETFREEUNIT()
    WRITE(testFile,'(A,I4.4,A)')"tests_",myRank,".out"
    OPEN(UNIT     = testUnit    ,&
         FILE     = testfile    ,&
         STATUS   = 'REPLACE'   ,&
         ACCESS   = 'SEQUENTIAL' )
  END IF

  !initialization phase
  dorestart=.FALSE.
  IF(PRESENT(RestartFile_in)) THEN
    dorestart=(LEN(TRIM(RestartFile_in)).GT.0)
  END IF
  IF(dorestart) CALL InitRestart(RestartFile_in)

  CALL InitOutput()
  CALL InitAnalyze()

  which_functional=GETINT('which_functional', Proposal=1 )
  CALL InitFunctional(functional,which_functional)


  CALL IgnoredStrings()

  CALL functional%InitSolution()
  StartTime=GetTime()
  SWRITE(Unit_stdOut,'(A,F8.2,A)') ' INITIALIZATION FINISHED! [',StartTime-StartTimeTotal,' sec ]'
  SWRITE(Unit_stdOut,fmt_sep)

  CALL functional%minimize()
  EndTime=GetTime()
  SWRITE(Unit_stdOut,'(A,2(F8.2,A))') ' FUNCTIONAL MINIMISATION FINISHED! [',EndTime-StartTime,' sec ], corresponding to [', &
       (EndTime-StartTime)/REAL(MaxIter,wp)*1.e3_wp,' msec/iteration ]'

  CALL FinalizeFunctional(functional)

  CALL FinalizeAnalyze()
  CALL FinalizeOutput()
  IF(dorestart) CALL FinalizeRestart()
  CALL FinalizeReadIn()


  ! do something
  IF(testlevel.GT.0)THEN
    SWRITE(UNIT_stdout,*)
    SWRITE(UNIT_stdOut,'(A)')"** TESTESTESTESTESTESTESTESTESTESTESTESTESTESTESTEST **"
    SWRITE(UNIT_stdout,*)
    n_warnings_occured=n_warnings_occured +nFailedMsg
    IF(nFailedMsg.GT.0)THEN
      SWRITE(UNIT_stdOut,'(A)')"!!!!!!!   SOME TEST(S) FAILED, see tests.out !!!!!!!!!!!!!"
    ELSE
      SWRITE(UNIT_stdOut,'(A)')"   ...   ALL IMPLEMENTED TESTS SUCCESSFULL ..."
    END IF !nFailedMsg
    SWRITE(UNIT_stdout,*)
    SWRITE(UNIT_stdOut,'(A)')"** TESTESTESTESTESTESTESTESTESTESTESTESTESTESTESTEST **"
    SWRITE(UNIT_stdout,*)
    CLOSE(testUnit)
  END IF !testlevel
  EndTimeTotal=GetTime()
  SWRITE(Unit_stdOut,fmt_sep)
  CALL DATE_AND_TIME(values=TimeArray) ! get System time
  SWRITE(UNIT_stdOut,'(A,I4.2,"-",I2.2,"-",I2.2,1X,I2.2,":",I2.2,":",I2.2)') &
    '%%% Sys date : ',timeArray(1:3),timeArray(5:7)
  SWRITE(Unit_stdOut,fmt_sep)
  IF(n_warnings_occured.EQ.0)THEN
    SWRITE(Unit_stdOut,'(A,F8.2,A)') ' GVEC SUCESSFULLY FINISHED! [',EndTimeTotal-StartTimeTotal,' sec ]'
  ELSE
    SWRITE(Unit_stdOut,'(A,F8.2,A,I8,A)') ' GVEC FINISHED! [',EndTimeTotal-StartTimeTotal,' sec ], WITH ' , n_warnings_occured , ' WARNINGS!!!!'
  END IF
  SWRITE(Unit_stdOut,fmt_sep)
  __PERFOFF('main')
  IF(MPIRoot) THEN
  __PERFOUT('main')
  END IF
END SUBROUTINE rungvec
END MODULE MODgvec_rungvec