Skip to content
Snippets Groups Projects
module_definitions.f90 16.1 KiB
Newer Older
  • Learn to ignore specific revisions
  • module definitions
    
    !=====[EDGE]====================================================================
    !this type is to store edges in a trianglulation
    ! it is used to update (in a generalized Delaunay sense)
    ! the triangulation of the 3D points on the surfaces
    ! for each edge:
    ! n1, n2 are the node numbers defining the edge
    ! t1, t2 are the triangle numbers on either side of the edge
    ! going from n1 to n2, t1 is to the left and t2 is to the right
    ! m1, m2 are the node numbers of the two other nodes making t1 and t2
    
          type edge
          integer n1,n2,m1,m2,t1,t2
          end type edge
    
    !=====[SHEET]===================================================================
    ! this type to store surfaces tracked by particles
    ! nsurface is number of particles
    ! x,y,z and xn,yn,zn are the coordinates and normals of the particles
    ! r,s are the 2D positions of the particles
    ! levelt is the resolution level (power of 2)
    ! itype is the type of surface (0 not foldable, 1 foldable)
    ! material is the material number
    
          type sheet
          integer nsurface,nt
          double precision,dimension(:),pointer::r,s,x,y,z,xn,yn,zn,u,v,w
    
          integer levelt,itype,material_main,material_spinup,material,surface_type
          integer closed
    
          double precision :: sp01,sp02,sp03,sp04,sp05,sp06,sp07,sp08,sp09,sp10,sp11
    
          double precision :: sp12,sp13,sp14,activation_time
    
          integer,dimension(:,:),pointer::icon
    
          logical rand,fixed_surf_spinup,fixed_surf
    
          integer leveloct
          double precision :: stretch,anglemax,anglemaxoctree
          integer :: criterion
          integer spread_surface_points
          end type sheet
    
    !=====[OCTREEV]=================================================================
    ! octreev type to store velocity/temperature octree
    ! octree is octree
    ! noctree is onctree max size
    ! nnode is number of nodes on the octree
    ! nleaves is number of leaves in the octree
    ! x,y,z are the coordinates of the nodes
    ! unode,vnode,wnode are the velocity conponents at the nodes
    
    ! wiso is the velocity z-component due to isostasy
    
    ! temp is temperature at the nodes
    ! icon is connectivity array between nodes
    
          type octreev
          integer,dimension(:),pointer::octree
          integer noctree,nnode,nleaves
          double precision,dimension(:),pointer::x,y,z,temp
    
          double precision,dimension(:),pointer::unode,vnode,wnode,wnodeiso
    
          double precision,dimension(:),pointer::pressure,spressure
    
          double precision,dimension(:),pointer::temporary_nodal_pressure
          integer,dimension(:,:),pointer::icon
          logical, dimension(:),pointer :: whole_leaf_in_fluid
          end type octreev
    
    !=====[OCTREELSF]===============================================================
    ! octreelsf type to store level set function
    ! octree is octree
    ! noctree is onctree max size
    ! nnode is number of nodes on the octree
    ! nleaves is number of leaves in the octree
    ! x,y,z are the coordinates of the nodes
    ! lsf is level set fuction values at the nodes
    ! icon is connectivity array between nodes
    
          type octreelsf
          integer,dimension(:),pointer::octree
          integer noctree,nnode,nleaves
          double precision,dimension(:),pointer::x,y,z,lsf
          integer,dimension(:,:),pointer::icon
          end type octreelsf
    
    !=====[OCTREESOLVE]=============================================================
    ! same as octreelsf except that this type can store nlsf lsfs
    ! total accuulated strain interpolated from a 3D cloud
    ! nodal velocities u,v,w, interpolated from octreev
    ! temperature velocities temp
    ! and bad faces
    
          type octreesolve
          integer,dimension(:),pointer::octree,kfix,kfixt
          integer noctree,nnode,nleaves,nlsf,nface
          double precision,dimension(:),pointer::x,y,z
          double precision,dimension(:),pointer::strain,pressure,spressure
          double precision,dimension(:),pointer::crit
    
          double precision,dimension(:),pointer::e2d,e2dp
    
          double precision,dimension(:),pointer::e3d
          double precision,dimension(:),pointer::lode
    
          double precision,dimension(:),pointer::u,v,w,wiso,temp
    
          double precision,dimension(:,:),pointer::lsf
    
          double precision,dimension(:),pointer :: eviscosity,q,dilatr,yield_ratio
          double precision,dimension(:),pointer :: frict_angle
    
          integer,dimension(:,:),pointer::icon,iface
          logical,dimension(:),pointer::is_plastic
    
          integer,dimension(:),pointer :: matnum
    
          end type octreesolve
    
    !=====[MATERIAL]================================================================
    ! type material to define different materials
    ! density is product of density by g
    ! viscosity is viscosity
    ! expon is the stress exponent in the expression for the viscosity
    ! activationenergy is the activation energy
    ! expansion is the coefficient of thermal expansion (note that the temperature is first multiplied by the
    ! temperature scale before it is multiplied by the expansion coefficient)
    ! penalty is incompressibility penalty factor
    ! plasticity_type is the type of plasticity (vM for von Mises)
    ! plasticity_parameters are the nine plasticity parameters
    
    ! mattrans are the xmin,xmax,ymin,ymax,zmin and zmax for material transitions
    ! transnum are the corresponding material numbers for the transition
    
    ! fviscosity is a scaling factor for power-law viscous behavior
    
    
          type material
          character(len=8) plasticity_type
    
          character(len=16) plasticity_ss_type_coh,plasticity_ss_type_phi
    
          character(len=16) fvisc_weak_type
    
          double precision :: plasticity_parameters(9),mattrans(6)
    
          double precision :: density,viscosity,penalty,expon,activationenergy
    
          double precision :: diffusivity,heat,expansion,fviscosity
    
          double precision :: fvisc_weak_onset,fvisc_weak_end,fvisc_weak_final
    
          end type material
    
    !=====[VOID]====================================================================
    ! type void to store information on node, leaves and faces that are in the void
    ! node=1 for nodes that are conpletely in the void (they are taken out of the equation set)
    ! leaf=1 for leaves that are completely in the void
    ! face=1 for faces that are completely in the void
    ! nnode is the number of nodes to be solved for (number of nodes not in the void)
    ! nleaves is the number of active leaves
    ! nface is the number of active faces
    ! rtf (restricted to full) is an array that provides the connectivity
    !   between the restricted and full
    !   node numbers (j=rtf(i) where i is restricted node number (1 to vo%nnode)
    !   and j is full node number (1 to nnode))
    ! ftr (full to restricted) is the complement
    ! influid is true for nodes that are in the fluid
    
          type void
          integer,dimension(:),pointer::node,leaf,face,ftr,rtf
          logical,dimension(:),pointer::influid
          integer nnode,nleaves,nface
          end type void
    
    !=====[CLOUD]===================================================================
    ! type cloud is to store volumetric cloud of points
    ! it is very simple and made of
    ! np : number of particles!
    ! x,y,z : location
    ! x0, y0, z0: original particle position (used to calculate F: deformation gradient
    !                                         and from it the total strain)
    ! strain: accumulated strain
    
    ! e2dp: effective strain rate from final iteration of previous time step
    
          integer np,ntag,matnum
    
          double precision,dimension(:),pointer :: x,y,z,x0,y0,z0,strain,lsf0,temp
          double precision,dimension(:),pointer :: press,e2dp
    
          integer,dimension(:),pointer::tag
          end type cloud
    
    !=====[TOPOLOGY]================================================================
    ! type topology is used to store the matrix topology for the solver
    
          type topology
          integer nheight,nheightmax
          integer,dimension(:),pointer::icol
          end type topology
    
    !=====[BOX]=====================================================================
    ! type box is used to store the geometrical information about a general box
    ! that is used to refine an octree in an arbitrary manner
    
          type box
          double precision x0,y0,z0,x1,y1,z1
          integer level
          end type box
    
    !=====[FACE]====================================================================
    ! type face is used to store the refinement informations on a given face of the cube
    
          type face
          integer level
          double precision t,b,l,r
          end type
    
    !=====[CROSS_SECTION]===========================================================
    ! this type is used to store all the parameters related to the cross sections 
    ! the 'flag' parameters toggle on/off the associated ield cross section
    ! scale is the global scale of the output postscript
    ! xyz =1,2,3 corresponding to a cross section at constant x,y, or z (resp.)
    
       type cross_section
          integer xyz
          double precision slice,scale
          logical flag_spress,flag_press,flag_e2d,flag_e3d,flag_crit,flag_grid,flag_lsf
          logical flag_u,flag_v,flag_w,flag_uvw,flag_q,flag_vfield,flag_strain
          logical flag_colour,flag_lode,flag_mu,flag_plastic,flag_velvect,flag_temp
          character(len=3) colormap
          integer ncolours
       end type
    
    !=====[PARAMETERS]==============================================================
    ! this contains all sorts of parameters used in the code, such as flags, 
    ! convergence parameters, octree levels, ... pretty much all parameters
    ! that are read in the input file.
    
    ! isobc is a boolean value used to turn on/off the use of the modified isostasy 
    ! boundary conditions (dwhipp 12/09)
    
    
       type parameters 
          logical compute_qpgram
          logical normaladvect
          logical excl_vol
          logical doDoRuRe
          logical ref_on_faces
          logical visualise_matrix
          logical renumber_nodes
          logical erosion
          logical adaptive_tol
          logical calculate_temp
          logical isostasy
          logical flexure
    
          logical bulkvisc
          logical init_e2d
    
          logical invariants_2d
    
          logical sedimentation
    
          logical remove_surf_pts
    
          logical sstemp_spinup
          logical in_spinup
          logical first_step
          logical first_iter
          logical first_nliter
          logical all_surf_fixed_spinup
          logical fixed_cloud_spinup
    
          logical materials_on_cloud
    
          integer nmpi
    
          integer levelcut
          integer levelapprox
          integer nstep
    
          integer noctreemax
          integer levelmax_oct
          integer leveluniform_oct
          integer smoothing_type
          integer debug
          integer npmin,npmax
          integer irestart
          integer refine_criterion
          integer niter_move
          integer griditer
    
          integer griditer_main
          integer griditer_spinup
    
          integer nonlinear_iterations
    
          integer nonlinear_iterations_main
          integer nonlinear_iterations_spinup
    
          integer mpe
          integer nmat
          integer nboxes
          integer nb_iter_nl_min
          integer nsections
          integer ns
          integer initial_refine_level
    
          integer baselevelx0
          integer baselevelx1
          integer baselevely0
          integer baselevely1
    
          integer sedimentation_type
    
          integer vermajor
          integer verminor
          integer verstat
          integer verrev
    
          integer svnrev
    
          integer,dimension(:),pointer::materialn
    
          double precision ztemp
          double precision dt
    
          double precision dt_main
          double precision dt_spinup
    
          double precision refine_ratio
          double precision courant
          double precision octree_refine_ratio
          double precision distance_exponent
          double precision tol
    
          double precision tol_main
          double precision tol_spinup
    
          double precision penalty
          double precision tempscale
          double precision zerosion
          double precision viscositymin,viscositymax
          double precision velocity_scale,length_scale
          double precision diffusion_erosion
          double precision fluvial_erosion
          double precision elastic_plate_thickness
          double precision density_difference
    
          double precision e2d0
    
          double precision plastic_stress_min
    
          double precision er_start
          double precision er_end
          double precision sed_start
          double precision sed_end
          double precision zaggrade_init
          double precision x_agg_start
          double precision x_agg_end
          double precision y_agg_start
          double precision y_agg_end
          double precision x_agg_sinus_amp
          double precision x_agg_sinus_wavelth
          double precision y_agg_sinus_amp
          double precision y_agg_sinus_wavelth
          double precision aggrade_rate
          double precision x_prog_start
          double precision x_prog_end
          double precision x_prog_length
          double precision z_prog_init
          double precision z_prog_fin
          double precision prog_rate_u
    
          double precision sstemp_viscosity_spinup
    
          character*128 restartfile
          character*40  infile
       end type
    
    !=====[INTERFACE SCANFILE]======================================================
    ! following is a general interface to read stuff from a file
    
          interface scanfile
            subroutine iscanfile (fnme,text,res,ires)
            character*(*) fnme,text
            integer,intent(out)::res
            integer,intent(out)::ires
            end subroutine iscanfile
            subroutine dscanfile (fnme,text,res,ires)
            character*(*) fnme,text
            double precision,intent(out)::res
            integer,intent(out)::ires
            end subroutine dscanfile
            subroutine cscanfile (fnme,text,res,ires)
            character*(*) fnme,text
            character*(*),intent(out)::res
            integer,intent(out)::ires
            end subroutine cscanfile
          end interface
    
    !=====[INTERFACE QSORT]=========================================================
    ! following is a general interface to sort an array of numbers
    
          interface qsort
    
            subroutine iqsort (array,n,perm)
            integer,intent(in) :: n
    
            integer,dimension(n),intent(inout)::array
            integer,dimension(n),intent(inout)::perm
    
            end subroutine iqsort
    
            subroutine rqsort (array,n,perm)
            integer,intent(in) :: n
    
            real,dimension(n),intent(inout)::array
            integer,dimension(n),intent(inout)::perm
    
            end subroutine rqsort
    
            subroutine dpqsort (array,n,perm)
            integer,intent(in) :: n
    
            real*8,dimension(n),intent(inout)::array
            integer,dimension(n),intent(inout)::perm
    
            end subroutine dpqsort
    
            subroutine iqsort_s (array,n)
            integer,intent(in) :: n
    
            integer,dimension(n),intent(inout)::array
    
            end subroutine iqsort_s
    
            subroutine rqsort_s (array,n)
            integer,intent(in) :: n
    
            end subroutine rqsort_s
    
    !        subroutine dpqsort_s (array,n)
    !        integer,intent(in) :: n
    !        real*8,dimension(n),intent(inout)::array
    !        end subroutine dpqsort_s
    
    
          end interface
    
    
    !=====[BC_DEFINITION]===========================================================
    ! type bc_definition is used to store values related to the applied boundary
    ! conditions, including the velocity values and certain geometric parameters
    
          type bc_definition
            character(len=3) :: bcorder
    
            character(len=40) :: bctype
    
            double precision,dimension(:,:),pointer :: zisodisp,zisoinc
    
            double precision :: bc_parameters(12),utrans,vtrans
            double precision :: ux0,ux1,vx0,vx1,wx0,wx1
            double precision :: uy0,uy1,vy0,vy1,wy0,wy1
    
            double precision :: uz0,uz1,vz0,vz1,wz0,wz1
            logical :: fixux0,fixux1,fixvx0,fixvx1,fixwx0,fixwx1
            logical :: fixuy0,fixuy1,fixvy0,fixvy1,fixwy0,fixwy1
            logical :: fixuz0,fixuz1,fixvz0,fixvz1,fixwz0,fixwz1
          end type bc_definition
    
    !=====[NEST_INFO]===============================================================
    ! type nest_info is used to store information related to the LS model and
    ! embedded SS nest
    !
    
    ! lsoutfile: file that is read to define the nest velocity B/Cs
    
    ! sselem[x,y,z]: ratio of ss element size to ls element size (should be 0-1)
    ! [x,y,z]minls: location of (0,0,0) in the SS model within the LS model
    
          type nest_info
    
            character(len=5)   :: ssoutdir
    
            character(len=128) :: lsoutfile
    
            double precision   :: sselemx,sselemy,sselemz,xminls,yminls,zminls
    
            double precision   :: sselemx0,sselemy0,sselemz0
    
    end module definitions