! ! File OrbScore2.f90 has these contents (in this order): ! ! BLOCK DATA BD1 ! BLOCK DATA BD2 ! PROGRAM OrbScore2 (which CONTAINS all its subprograms, although any global references are NOT intended) ! !================================================================================== BLOCK DATA BD1 ! Define "weight" (Gaussian integration weights) of the ! seven integration points in each element, defined by internal ! coordinates "points(3, 7)", where points(1:3, m) = s1:s3 OF ! integration point number m. ! Because all of these arrays are functions of internal ! coordinates, they are not affected by scaling or deformation of ! the elements. IMPLICIT NONE DOUBLE PRECISION points, weight COMMON / S1S2S3 / points COMMON / WgtVec / weight DIMENSION points(3, 7), weight(7) ! "points" contains the internal coordinates (s1, s2, s3) of the 7 ! Gaussian integration points (for area integrals) of the ! triangular elements. "points" is also the set of nodal functions ! for unprojected scalar quantities within an element. DATA points / & & 0.3333333333333333D0, 0.3333333333333333D0, 0.3333333333333333D0, & & 0.0597158733333333D0, 0.4701420633333333D0, 0.4701420633333333D0, & & 0.4701420633333333D0, 0.0597158733333333D0, 0.4701420633333333D0, & & 0.4701420633333333D0, 0.4701420633333333D0, 0.0597158733333333D0, & & 0.7974269866666667D0, 0.1012865066666667D0, 0.1012865066666667D0, & & 0.1012865066666667D0, 0.7974269866666667D0, 0.1012865066666667D0, & & 0.1012865066666667D0, 0.1012865066666667D0, 0.7974269866666667D0 / ! "weight" is the Gaussian weight (for area integrals) of the 7 ! integration points in each triangular element. DATA weight / 0.2250000000000000D0, & & 0.1323941500000000D0, 0.1323941500000000D0, 0.1323941500000000D0, & & 0.1259391833333333D0, 0.1259391833333333D0, 0.1259391833333333D0 / END BLOCK DATA BD1 BLOCK DATA BD2 ! Define "fPhi" (nodal functions) and "fGauss" (Gaussian integration ! weights) at the 7 integration points in each fault element, ! defined by internal coordinate "fPoint(m = 1:7)" ! which contains the relative positions ! (fractional length) of the integration points. ! Because all of these arrays are functions of internal ! coordinates, they are not affected by scaling or deformation of ! the elements. IMPLICIT NONE DOUBLE PRECISION fPhi, fPoint, fGauss COMMON / SFault / fPoint COMMON / fPhis / fPhi COMMON / FGList / fGauss DIMENSION fPhi(4, 7), fPoint(7), fGauss(7) ! "fPoint" contains the seven integration point locations for the fault ! elements. Each value gives a position as a fraction of total length ! measured from node1 to node2 (of array "nodeF"). DATA fPoint / & & 0.0254461D0, & & 0.1292344D0, & & 0.2970774D0, & & 0.5000000D0, & & 0.7029226D0, & & 0.8707656D0, & & 0.9745539D0 / ! "fGauss" contains the seven corresponding weight factors for use in ! line integrals. DATA fGauss / & & 0.0647425D0, & & 0.1398527D0, & & 0.1909150D0, & & 0.2089796D0, & & 0.1909150D0, & & 0.1398527D0, & & 0.0647425D0 / ! "fPhi" contains the values of the 4 nodal functions (one per node) ! at each of these 7 integration points in the fault element. ! A special convention is that the nodel function of node 3 ! is the negative of that for node 2, while the nodal function ! for node 4 is the negative of that for node 1. This simplifies ! many expressions in which we would otherwise have to have ! a separate factor of +1.0D0 or -1.0D0 for the two sides of the fault. DATA fPhi / & & 0.9745539D0, 0.0254461D0, -0.0254461D0, -0.9745539D0, & & 0.8707656D0, 0.1292344D0, -0.1292344D0, -0.8707656D0, & & 0.7029226D0, 0.2970774D0, -0.2970774D0, -0.7029226D0, & & 0.5000000D0, 0.5000000D0, -0.5000000D0, -0.5000000D0, & & 0.2970774D0, 0.7029226D0, -0.7029226D0, -0.2970774D0, & & 0.1292344D0, 0.8707656D0, -0.8707656D0, -0.1292344D0, & & 0.0254461D0, 0.9745539D0, -0.9745539D0, -0.0254461D0 / END BLOCK DATA BD2 !================================================================================== PROGRAM OrbScore2 ! Compares output from -Shells- with data from: ! -geodetic networks, ! -stress directions, ! -fault slip rates, ! -seafloor spreading rates, ! -seismicity, ! -upper-mantle seismic anisotropy, ! and reports both comparison tables, and summary scalar misfits. !=========== Part of the -Shells- package of programs =========== ! Given a finite-element grid file in the format produced by -OrbWin- ! (or -OrbWeave-) and renumbered by -OrbNumber-, with nodal data ! added by -OrbData5- (or -OrbData-), and node-velocity output from -Shells-, ! computes a variety of misfit measures for the realism of the results. ! NOTES: *Does not contain subprograms Viscos or Diamnd, hence independent ! of changes made in May 1998, and equally compatible ! with -Old_SHELLS- or with improved -Shells-. ! (However, subprogram Mohr IS used to find brittle/ductile ! transition depths in fault elements.) ! *Reads any optional Lithospheric Rheology index numbers ! (LR#s) in different fault and/or continuum elements ! (a feature of Shells_v5.0+). Uses these in CALL Mohr. ! by ! Peter Bird ! Department of Earth, Planetary, and Space Sciences, ! University of California, Los Angeles, California 90095-1567 ! (C) Copyright 1994, 1998, 1999, 2000, 2001, 2002, 2005, 2006, ! 2015, 2018, 2019, 2021 ! by Peter Bird and the Regents of ! the University of California. ! (For version date, search for "Version of" below) ! This program was developed with support from the University of ! California, the United States Geologic Survey, the National ! Science Foundation, and the National Aeronautics and Space ! Administration. ! It is freeware, and may be copied and used without charge. ! It may not be modified in a way which hides its origin ! or removes this message or the copyright message. ! It may not be resold for more than the cost of reproduction ! and transmission. ! !--------------------------------------------------------------------- ! *** ADDITIONAL SOFTWARE REQUIRED TO SOLVE A LINEAR SYSTEM: *** ! ! Using one routine (dsysv) from the LAPACK (Linear Analysis Package) ! portion of the Intel Math Kernel Library (MKL): USE MKL95_LAPACK ! in file LAPACK.f90, available from Peter Bird (I just renamed LAPACK.inc). ! ! {Note that program lines associated with each of these solvers are just ! commented-out, not removed, when the solver package is changed.} !--------------------------------------------------------------------- ! *** ADDITIONAL Fortran 90 code needed to compile: USE DSphere ! provided by Peter Bird as file DSphere.f90 USE DDislocation ! provided by Peter Bird as file DDislocation.f90 USE DIcosahedron ! provided by Peter Bird as file DIcosahedron.f90 USE DWeighting ! provided by Peter Bird as file DWeighting.f90 !--------------------------------------------------------------------- ! ! DIFFERENCES between OrbScore and OrbScore2: ! ! (1) OrbScore was in fixed-format, due to its FORTRAN 77 heritage. But, ! OrbScore2 has been (nominally, minimally) converted to free-format ! Fortran 90 code, which is clearer and easier to read and modify. ! (2) OrbScore2 uses MODULE DWeighting to weight geodesy, stress, ! & anisotropy by geographic area (but not spreading-rates ! or fault-slip0rates or seismicity): ! N.B. MODULE DWeighting has the PRIVATE attribute to shield all ! of its internal variables and arrays from accidental ! direct access. This is critical because MODULE DWeighting ! creates and stores a *different* F-E grid (a uniform subdivision ! of the icosahedron) for purposes of computing data weights. ! (3) IF (pltGeo) a new output file provides errors in predicted velocities ! of geodetic benchmarks, for plotting (e.g., with -FiniteMap-). ! IF (pltCos) a new output file provides coseismic part of model velocities ! of geodetic benchmarks, for plotting (e.g., with -FiniteMap-). ! IF (pltInt) a new output file provides interseismic part of model velocities ! of geodetic benchmarks, for plotting (e.g., with -FiniteMap-). ! (4) OrbScore2 includes scoring for match between upper-mantle seismic ! anisotropy (specifically, fast-direction azimuths from SKS splitting) ! and basal shear traction vectors computed from a torque report, ! in runs with iConve = 6. This scoring is area- and delay-weighted. ! N.B. This code, added 2006.11 for the Earth5 project, still needs ! to be generalized to other iConve options for basal shear traction! ! (5) OrbScore2 reads and uses any (optional) Lithospheric Rheology index ! codes following the definitions of continuum-triangle and/or ! linear-fault elements (a Shells_v5.0+ feature); OrbScore does not. !--------------------------------------------------------------------- IMPLICIT NONE ! All variable names must be declared & typed. !--------------------------------------------------------------------- ! PARAMETER (ARRAY-SIZE) STATEMENTS ! Set the following parameters at least as large as your problem: ! nPlate = Number of pLates in PB2002 model of Bird [2003]: INTEGER, PARAMETER :: nPlate = 52 ! nPBnd = Maximim number of boundary points (per plate) in PB2002 model of Bird [2003]: INTEGER, PARAMETER :: nPBnd = 1250 ! maxNod = Maximun number of nodes (includes both "real" & "fake"): INTEGER, PARAMETER :: maxNod = 17000 ! maxEl = Maximun number of continuum elements (triangles): INTEGER, PARAMETER :: maxEl = 27000 ! maxFEl = Maximun number of fault elements (lines); INTEGER, PARAMETER :: maxFEl = 3000 ! maxBN = Maximum number of boundary nodes: INTEGER, PARAMETER :: maxBN = 2000 ! maxAtP = Maximum number of nodes which may overlap at a fault- ! intersection point: INTEGER, PARAMETER :: maxAtP = 20 ! maxGeo = Maximum number of geodetic benchmarks at which geodetic ! horizontal velocity vectors are provided (for scoring): INTEGER, PARAMETER :: maxGeo = 20000 ! maxStr = Maximum number of points at which most-compressive ! horizontal principal stress direction is provided ! (for scoring): INTEGER, PARAMETER :: maxStr = 10000 ! maxSKS = maximum number of points at which fast-polarization ! ("phi") azimuth and splitting time of SKS waves are ! provided (for scoring). INTEGER, PARAMETER :: maxSKS = 10000 ! maxDat = Maximum number of sites with data on fault slip rate: INTEGER, PARAMETER :: maxDat = 5000 ! maxMOR = Maximum number of points at which seafloor-spreading ! rates are provided (for scoring): INTEGER, PARAMETER :: maxMOR = 800 ! maxEqs = Maximum number of earthquakes in catalog used for scoring: INTEGER, PARAMETER :: maxEqs = 30000 ! maxAdj = Maximum number of (longitude/latitude/velocity-datum/ ! predicted-velocity) points passed to -Adjust- for application ! of a rigid-body rotation-rate that minimizes the weighted ! RMS velocity difference: INTEGER, PARAMETER :: maxAdj = 188000 ! "piO180" is Pi/180 (conversion from degrees to radians): REAL*8, PARAMETER :: piO180 = 0.0174532925199433D0 !--------------------------------------------------------------------- ! TYPE & DIMENSION STATEMENTS associated with supporting Lithospheric Rheology index #s: CHARACTER*80 :: longer_line, shorter_line INTEGER :: LRi, LRn ! <== LRn = highest LRi# found in the .FEG file; used in allocating arrays just below: LOGICAL, DIMENSION(:), ALLOCATABLE :: LR_is_defined, LR_is_used ! (0:LRn) REAL*8 :: d_fFric, d_cFric, d_Biot, d_Byerly, d_eCreep ! <=== scalars in default rheology REAL*8 :: d_aCreep(2), d_bCreep(2), d_cCreep(2), d_dCreep(2) ! <=== default rheologic constants for creep which have different crust:mantle values REAL*8, DIMENSION(:), ALLOCATABLE :: LR_set_fFric, LR_set_cFric, LR_set_Biot, LR_set_Byerly, & ! (0:LRn) & LR_set_eCreep REAL*8, DIMENSION(:,:), ALLOCATABLE :: LR_set_aCreep, LR_set_bCreep, LR_set_cCreep, LR_set_dCreep ! (1:2, 0:LRn) !--------------------------------------------------------------------- ! DIMENSION STATMENTS ! DIMENSIONs using PARAMETER maxNod: LOGICAL :: checkN REAL*8 :: atnode, cooling_curvature, density_anomaly, dQdTdA, eDotNC, eDotNM, elev, & & tLNode, uvecN, xNode, yNode, zMNode DIMENSION atnode (maxNod), checkN (maxNod), & & cooling_curvature(maxNod), density_anomaly(maxNod), & & dQdTdA (maxNod), & & eDotNC (maxNod), eDotNM (maxNod), & & elev (maxNod), & & tLNode (maxNod), & & uvecN (3, maxNod), v (2, maxNod), & & xNode (maxNod), yNode (maxNod), zMNode (maxNod) ! DIMENSIONs using PARAMETER maxEl: INTEGER :: continuum_LRi(maxEl), nodes LOGICAL :: checkE REAL*8 :: area, CartR, detJ, dxs, dys, dxsp, dysp, & & eDotEC, eDotEM, eRate, fpsfer, sita, tLInt, xIP, yIP, zMoho DIMENSION area (maxEl), & & CartR (3, 7, maxEl), checkE (maxEl), detJ (7, maxEl), & & dxs (2, 2, 3, 7, maxEl), dys (2, 2, 3, 7, maxEl), & & dxsp (3, 7, maxEl), dysp (3, 7, maxEl), edgeTs (3, maxEl), & & eDotEC (maxEl), eDotEM (maxEl), & & eRate (3, 7, maxEl), & & fpsfer (2, 2, 3, 7, maxEl), & & nodes (3, maxEl), & & sita (7, maxEl), & & tLInt (7, maxEl), & & xIP (7, maxEl), yIP (7, maxEl), & & zMoho (7, maxEl) ! DIMENSIONs using PARAMETER maxFEl: INTEGER :: fault_LRi(maxFEl), nodeF LOGICAL :: checkF REAL*8 :: fc, fDip, fimudz, fLen, fpeaks, & & fPFlt, fArg, ftstar, offset, zTranF DIMENSION checkF (maxFEl), edgeFS (2, maxFEl), & & fc (2, 2, 7, maxFEl), fDip (2, maxFEl), & & fimudz (7, maxFEl), fLen (maxFEl), fpeaks (2, maxFEl), & & fPFlt (2, 2, 2, 7, maxFEl), fSlips (maxFEl), & & fArg (2, maxFEl), ftstar (2, 7, maxFEl), nodeF (4, maxFEl), & & offset (maxFEl), zTranF (2, maxFEl) ! DIMENSIONs using PARAMETER maxBN: INTEGER :: nodCon DIMENSION nodCon(maxBN) ! DIMENSIONs of (2) refer to crust/mantle-lithosphere: REAL*8 :: alphaT, conduc, & & radio, rhoBar, temLim DIMENSION alphaT(2), conduc(2), & & radio(2), rhoBar(2), temLim(2) ! DIMENSIONs using PARAMETER maxAtP: INTEGER :: list DIMENSION list (maxAtP) ! DIMENSIONs using PARAMETER nPlate: REAL*8 :: omega DIMENSION names(nPlate), omega(3, nPlate) INTEGER, DIMENSION(nPlate) :: nBoundaryPoints REAL*8, DIMENSION(nPlate, nPBnd) :: pLat, pLon ! DIMENSIONs using PARAMETER maxGeo: CHARACTER*20 :: geoTag REAL*8 :: geoPhi, geoThe, geoVel, geoSig, geoAzi, & & geoUTh, geoUPh DIMENSION geoTag(maxGeo), geoPhi(maxGeo), geoThe(maxGeo), & & geoVel(maxGeo), geoSig(maxGeo), geoAzi(maxGeo), & & geoUTh(maxGeo), geoUPh(maxGeo) ! DIMENSIONs using PARAMETER maxStr: CHARACTER*5 :: strTag(maxStr) REAL*8 :: strPhi, strThe, strQua, strArg DIMENSION strPhi(maxStr), strThe(maxStr), & & strQua(maxStr), strArg(maxStr), strReg(maxStr) ! DIMENSIONs using PARAMETER maxSKS: CHARACTER*5, DIMENSION(maxSKS) :: SKS_tag ! 5-byte ID code for datum REAL*8, DIMENSION(maxSKS) :: SKS_theta ! colatitude, in radians REAL*8, DIMENSION(maxSKS) :: SKS_phi ! colatitude, in radians REAL*8, DIMENSION(maxSKS) :: SKS_argument ! phi or fast direction, in radians counterclockwise from S REAL*8, DIMENSION(maxSKS) :: SKS_delay ! SKS splitting time, in s INTEGER, DIMENSION(maxSKS) :: SKS_iPlate ! integer identifying PB2002 plate which contains this datum REAL*8, DIMENSION(2, maxSKS) :: basal_shear_tractions ! (theta = S, phi = E) components, in Pa ! DIMENSIONs using PARAMETER maxDat: REAL*8 :: delVZ, fltLat, fltLon, rLat DIMENSION delVZ(2, maxDat), fName(maxDat), & & fltLat(maxDat), fltLon(maxDat), & & rLat(2, maxDat) ! DIMENSIONs using PARAMETER maxMOR: CHARACTER*5 :: MORTag(maxMOR) REAL*8 :: MORPhi(maxMOR), MORThe(maxMOR), & & MORVel(maxMOR), MORSig(maxMOR) ! DIMENSIONs using PARAMETER maxEqs: REAL*8 :: eqELon, eqNLat, eqMag DIMENSION eqELon(maxEqs), eqNLat(maxEqs), eqMag(maxEqs) ! DIMENSIONs using PARAMETER maxAdj: REAL*8 :: data, predic, longterm_at_GPS, eLon ! N.B. nLat is already declared REAL*8 DIMENSION data(2, maxAdj), predic(2, maxAdj), longterm_at_GPS(2, maxAdj), & & eLon (maxAdj), nLat (maxAdj) ! DIMENSIONs of fixed size: REAL*8 :: CartVs, cUvec, equVec, rA, rB, rS, tauMax, tempV, & & uP, uT, vf, voff, vt DIMENSION CartVs(3, 3), cUvec(3), & & equVec(3), & & rA(3), rB(3), rS(3), tauMax(2), tempV(3), & & uP(3), uT(3), vf(3), voff(3), vt(3) ! Following statement to agree with BLOCK DATA BD1: DIMENSION points(3, 7), weight(7) ! Following statement to agree with BLOCK DATA BD2: DIMENSION fPhi(4, 7), fPoint(7), fGauss(7) ! ALLOCATABLE array(s) added for OrbScore2, to work with MODULE DWeighting: REAL*8, DIMENSION(:), ALLOCATABLE :: weights !--------------------------------------------------------------------- ! COMMON STATEMENTS ! Note: Un-named COMMON passes INTEGER variables used in the ! INTEGER FUNCTION "IndexK", to avoid passing these same ! through long sequences of subprograms. INTEGER :: lda, md COMMON lda, md ! Named COMMON BLOCKs BD1 & BD2 hold the positions, ! weights, and nodal function values at the integration points ! in the elements (triangular elements in BLOCK DATA BD1, and ! fault elements in BLOCK DATA BD2). ! Entries corresponding to BD1: COMMON / S1S2S3 / points COMMON / WgtVec / weight ! Entries corresponding to BD2: COMMON / SFault / fPoint COMMON / fPhis / fPhi COMMON / FGList / fGauss !-------------------------------------------------------------------- ! DATA STATEMENTS ! "floats" is a switch controlling whether subprogram Adjust ! will be used to minimize the kinetic energy (RMS velocity) ! of velocity-vector lists in two contexts: ! (1) In computing RMS velocities of nodes; ! (2) In computing mismatch between geodetic and model velocities. ! For normal use, set floats =.TRUE. ! If you are *SURE* that your geodetic data are expressed in the ! same velocity reference frame as is used to express velocities ! in your finite element model, select floats = .FALSE. for a ! more discriminating (less forgiving) test of model quality. LOGICAL, PARAMETER :: floats = .TRUE. ! OPTION for a detailed table of offset-rate misfits ! (possibly useful when running OrbScore2 interactively): LOGICAL, PARAMETER :: table_of_offset_rates = .TRUE. ! "oezOpi" is 180/Pi (conversion from radians to degrees): REAL*8, PARAMETER :: oezOpi = 57.2957795130823D0 ! "secPYr" is the number of seconds in a year: REAL*8, PARAMETER :: secPYr = 3.155760D7 ! "unitL" is the number of millimeters per unit length in the ! system of units in use on Fortran devices iUnitG (finite ! element grid), and device iUnitP (parameters), ! and iUnitV (velocity solutions). For example, if these are in ! SI units (meters), then unitL = 1000.0D0 REAL*8, PARAMETER :: unitL = 1000.0D0 ! "unitT" is the number of years per unit of time in the ! system of units in use on Fortran devices iUnitG (finite ! element grid), and devide iUnitP (parameters), ! and iUnitV (velocity solutions). For example, if these are in ! SI units (seconds), then unitT = 3.1688087D-08 REAL*8, PARAMETER :: unitT = 3.1688087D-08 ! "dipMax" is the maximum dip (from horizontal, in degrees) for a ! fault element to be treated as a dip-slip fault, with two degrees ! of freedom per node-pair. At steeper dips, the degree of freedom ! corresponding to opening or convergence of the opposite sides is ! eliminated by a constraint equation, and the fault is treated as ! a vertical strike-slip fault. This arbitrary limit is necessary ! because the equations for dip-slip faults become singular as the ! dip approaches 90 degrees. In practice, it is best to specify dips ! as either (1) vertical, or (2) clearly less than "dipMax", within ! each fault element. If the dip varies within an element in such a ! way that it passes through this limit within the element, then ! the representation of that fault element in the global equations ! may be inaccurate. REAL*8, PARAMETER :: dipMax = 75.0D0 ! "subDip" is the maximum dip (from horizontal, in degrees) for a ! fault in a whole-Earth model (sphere = .TRUE.) to be treated as ! a subduction zone (in which case, the footwall requires ! boundary conditions). ! In all models, faults with less than this dip have the down-dip ! integral of traction limited to "tauMax". tauMax is an array ! of two values, for oceanic and continental subduction zones, ! respectively. If such limits are not wanted, then ! the tauMax values can just be set to very large numbers. REAL*8, PARAMETER :: subDip = 19.0D0 ! THE FOLLOWING ARE THE FORTRAN INPUT AND OUTPUT DEVICE NUMBERS: ! "iUnitT"= device number associated with text output, including ! status and warning messages. ! (NOTE! On some systems, system error messages are always ! output on device 6. If so, then iUnitT should be 6.) INTEGER, PARAMETER :: iUnitT = 6 ! "iUnitG"= device number associated with the .FEG grid input file: INTEGER, PARAMETER :: iUnitG = 1 ! "iUnitP"= device number associated with the .IN parameter input file. INTEGER, PARAMETER :: iUnitP = 2 ! "iUnitV"= device number associated with the v_____.OUT velocity ! solution (at the nodes) file: INTEGER, PARAMETER :: iUnitV = 3 ! "iUnitB"= device number associated with the geodetic file of ! benchmark velocities: INTEGER, PARAMETER :: iUnitB = 11 ! "iUnitS"= device number associated with the file of maximum ! horizontal principal stress directions: INTEGER, PARAMETER :: iUnitS = 12 ! "iUnitD"= device number associated with the geological slip rate file ! (in .DIG format, with two title lines per point): INTEGER, PARAMETER :: iUnitD = 13 ! "iUnitM"= device number associated with the file of mid-ocean ! spreading rates: INTEGER, PARAMETER :: iUnitM = 14 ! "iUnitE"= device number associated with the file of earthquakes ! to be used for seismicity scoring: INTEGER, PARAMETER :: iUnitE = 15 ! "iUnitO"= device number to be used for outputting the ! .FEG file with scalar strain-rates at nodes ! based on the earthquake catalog from iUnitE: INTEGER, PARAMETER :: iUnitO = 21 ! "iUnitZ"= device number to be used for outputting the ! .FEG file with scalar strain-rates at nodes ! based on the model (after healing faults and smoothing): INTEGER, PARAMETER :: iUnitZ = 22 ! "iUnitX"= Device number to be used for outputting the ! smoothed (interseismic) velocities of each node ! in v_____.out format, for use with FiniteMap. ! (This special file is only created IF(EXPLOD), ! which is true only if the "benchmarks" were ! created by Explode3 and if they match one-for- ! one with nodes of the .feg file.) INTEGER, PARAMETER :: iUnitX = 23 ! "iUnitY"= Device number to be used for outputting the ! error in geodetic velocity predictions, in same ! .gps format as the input dataset. ! Note that this occurs only IF (pltGeo) .AND. ! (some geodetic data are provided for scoring). INTEGER, PARAMETER :: iUnitY = 31 LOGICAL, PARAMETER :: pltGeo = .TRUE. ! "iUnitC"= Device number to be used for outputting the model ! coseismic part of geodetic velocity predictions, ! in same .gps format as the input dataset. ! Note that this occurs only IF (pltCos) .AND. ! (some geodetic data are provided for scoring). INTEGER, PARAMETER :: iUnitC = 32 LOGICAL, PARAMETER :: pltCos = .TRUE. ! "iUnitI"= Device number to be used for outputting the model ! interseismic part of geodetic velocity predictions, ! in same .gps format as the input dataset. ! Note that this occurs only IF (pltInt) .AND. ! (some geodetic data are provided for scoring). INTEGER, PARAMETER :: iUnitI = 37 LOGICAL, PARAMETER :: pltInt = .TRUE. ! "iUnitK"= Device number to be used for reading file with ! upper-mantle anisotropy data, in the form of fast- ! polarization ("phi") azimuths from SKS splitting, ! together with SKS splitting times in s. INTEGER, PARAMETER :: iUnitK = 33 ! "iUnitQ"= Device number to be used for reading file with ! torque report (obtained by running Shells version ! of summer 2006 or later). This is used in scoring ! the anisotropy data read through iUnitK (above). INTEGER, PARAMETER :: iUnitQ = 34 ! "iUnitF"= Device number to be used for reading file with ! outlines of the plates (e.g., PB2002_plates.dig). ! This is used in scoring the anisotropy data read ! through iUnitK (above). INTEGER, PARAMETER :: iUnitF = 35 ! iUnitLR = device number associated with non-default Lithospheric Rheologies: INTEGER, PARAMETER :: iUnitLR = 36 ! PB2002 plate names [Bird, 2003, G**3, Table 1]: DATA names / & & 'AF', 'AM', 'AN', 'AP', 'AR', 'AS', 'AT', 'AU', 'BH', 'BR', & & 'BS', 'BU', 'CA', 'CL', 'CO', 'CR', 'EA', 'EU', 'FT', 'GP', & & 'IN', 'JF', 'JZ', 'KE', 'MA', 'MN', 'MO', 'MS', 'NA', 'NB', & & 'ND', 'NH', 'NI', 'NZ', 'OK', 'ON', 'PA', 'PM', 'PS', 'RI', & & 'SA', 'SB', 'SC', 'SL', 'SO', 'SS', 'SU', 'SW', 'TI', 'TO', & & 'WL', 'YA' / INTEGER, PARAMETER :: iPAfri = 1 ! Index number of Africa plate in this model. ! Following rotation vectors in Cartesian (x,y,z) components, ! with units of radians per million years: ! [Bird, 2003, G**3, Table 1] DATA ((omega(i, j), i = 1, 3), j = 1, nPlate) / & & 0.002401, -0.007939, 0.013892, & ! 1 & 0.000949, -0.008643, 0.013725, & ! 2 & 0.000689, -0.006541, 0.013676, & ! 3 & 0.002042, -0.013153, 0.008856, & ! 4 & 0.008570, -0.005607, 0.017497, & ! 5 & 0.000148, -0.003070, 0.010915, & ! 6 & 0.015696, 0.002467, 0.023809, & ! 7 & 0.009349, 0.000284, 0.016252, & ! 8 & 0.000184, 0.005157, 0.001150, & ! 9 & -0.000871, -0.002268, 0.002507, & ! 10 & -0.019124, 0.030087, 0.010227, & ! 11 & 0.011506, -0.044526, 0.007197, & ! 12 & 0.001688, -0.009048, 0.012815, & ! 13 & 0.003716, -0.003791, 0.000949, & ! 14 & -0.008915, -0.026445, 0.020895, & ! 15 & -0.061175, 0.005216, -0.013755, & ! 16 & 0.070136, 0.160534, 0.094328, & ! 17 & 0.000529, -0.007235, 0.013123, & ! 18 & -0.083251, -0.002464, -0.014923, & ! 19 & 0.016256, 0.089364, 0.015035, & ! 20 & 0.008180, -0.004800, 0.016760, & ! 21 & 0.006512, 0.003176, 0.005073, & ! 22 & 0.108013, 0.299461, 0.230528, & ! 23 & 0.033318, -0.001813, 0.036441, & ! 24 & -0.013835, 0.008245, 0.015432, & ! 25 & -0.777844, 0.440872, -0.047437, & ! 26 & 0.001521, 0.007739, 0.013437, & ! 27 & 0.038223, -0.058291, 0.013679, & ! 28 & 0.001768, -0.008439, 0.009817, & ! 29 & -0.004336, 0.003769, -0.000402, & ! 30 & 0.000111, -0.006361, 0.010449, & ! 31 & 0.044913, -0.009546, 0.010601, & ! 32 & -0.055342, -0.010890, 0.006794, & ! 33 & -0.000022, -0.013417, 0.019579, & ! 34 & 0.001041, -0.008305, 0.012143, & ! 35 & -0.026223, 0.020184, 0.037208, & ! 36 & 0.000000, 0.000000, 0.000000, & ! 37 & -0.000040, -0.009291, 0.012815, & ! 38 & 0.012165, -0.012510, -0.000366, & ! 39 & -0.019183, -0.070604, 0.036798, & ! 40 & 0.000472, -0.006355, 0.009100, & ! 41 & 0.121443, -0.078836, 0.027122, & ! 42 & 0.001117, -0.007434, 0.008534, & ! 43 & -0.000833, -0.006701, 0.013323, & ! 44 & 0.001287, -0.008754, 0.014603, & ! 45 & -0.017196, 0.017186, 0.008623, & ! 46 & 0.003201, -0.010440, 0.015854, & ! 47 & 0.023380, -0.019369, -0.010465, & ! 48 & -0.009400, 0.023063, 0.008831, & ! 49 & 0.142118, 0.005616, 0.078214, & ! 50 & -0.016831, 0.018478, 0.010166, & ! 51 & -0.000836, -0.006169, 0.016274 / ! 52 !--------------------------------------------------------------------- ! STATEMENT FUNCTION(S): ! Interpolation within one fault element: REAL*8 :: fhival, s, x1, x2 fhival(s, x1, x2) = x1 + s * (x2 - x1) !--------------------------------------------------------------------- ! NON-Array (scalar) specifications: ! CHARACTER*132 :: gpsFMT CHARACTER*100 :: fName CHARACTER*80 :: enctit, line, title1, title2, title3 CHARACTER*20 :: tag CHARACTER*15 :: frame CHARACTER*10 :: flag10 CHARACTER*8 :: c8r1, c8r2, c8rerr, c8v1, c8v2, c8verr CHARACTER*7 :: c7rcha, c7vcha CHARACTER*4 :: outcom CHARACTER*3 :: c3 CHARACTER*2 :: names, reg, strreg INTEGER :: i, iConve, iDepth, iEle, iHour, iMohr, ios, iPlate, iPVRef, iSecon, iTenth, & & j, j1, j2, jYear, k, kDay, kp1, & & m, maxItr, mIP, minute, month, mxDOF, mxNode, mxEl, mxFEl, mxBn, mxStar, & & mxGeo, mxStr, mxSKS, mxGSR, mxMOR, mxAdj, & & n, n1, n2, n3, n4, n_AB_data, nB, nBadRe, nBland, nCond, nFakeN, nFData, nFl, nGSDat, nInSum, nJ1, nJ2, nL1, nL2, nL3, nL4, & & node, node1, node2, nRealN, nSKS_bad, nSKS_used, numAdj, numDen, numEl, numEqs, numGeo, numMOR, numNod, numSKS, numStr, n1000, & & subdivision LOGICAL :: brief, everyP, e1Part, e2Part, eZPart, & & haveNV, log_strike_adjustments, showis, skipBC, sphere, vertic LOGICAL, DIMENSION(nPlate) :: slab_q ! Note: The following arrays could be compressed with LOGICAL*1: LOGICAL :: edgeTs, edgeFS, explod, fSlips ! N.B. DOUBLE PRECISION specifications predate the current REAL*8 ! specifications, although the meaning is the same. DOUBLE PRECISION :: s_dp, shrink, v ! Following statement to agree with BLOCK DATA BD1: DOUBLE PRECISION :: points, weight ! Following statement to agree with BLOCK DATA BD2: DOUBLE PRECISION :: fPhi, fPoint, fGauss REAL*8 :: a, aBadRe, aCos, angle, anisotropy, arc, arg, aRegime_denominator, aSum1, aSum1D, aSum2, aSum2D, aveC, aveM, azim, azimuth, & & b, bad, badMma, badSig, bAzim, bigB, & & chance, closes, constr, correl, cosA, cosB, croSiz, crossP, & & deltaT, delVX, delVY, dExtra, dip, dot, dVAP, dVAT, dVBP, dVBT, dVP, dVT, dX, dY, & & e1, e11h, e12h, e1h, e2, e22h, e2h, e2me1, e2meZ, e3, e3azim, eLarge, eLonDe, eLonP, & & enB, eqM0, eqMB, eqMS, eqMW, equPar, err, etaMax, exx, exy, eyy, eZ, eZme1, & & factor, fLong, fMuMax, frac, fSave, & & geodes, gMean, gradie, gVMma, & & mvmma, & & nLat, nLatDe, nLatP, & & offMax, OKDelV, OKToQt, oneKm, opens, & & pAzim, perBad, phi, pure_bad, probRL, probVZ, & & r, r2, r2flt, r2min, r2n1, r2n2, radius, rate, refStr, regime, relVZ, rErr, rho, rhoAst, rhoH2O, rLt, rMin, rmsv, & & s1, s2, s3, seismi, short, side, sigma, sigma2, simple_bad, SIMu, sinist, slide, slpErr, smallA, smallB, smallC, spread, sprMma, st, stress, & & sum, sum0s, sum1, sum1d, sum1mm, sum1s, sum2, sum2d, sum2mm, sum2s, sumB,sumB2, sumN, sumNmm, sumNs, sumSid, & & t, t1, t2, tAdiab, test, theLat, theLon, theta, thick, time1, time2, toler, trHMax, tSec, tSurf, & & u1x, u1y, u2x, u2y, unitBX, unitBY, unitX, unitY, & & varC, varM, vEMmpa, vESigm, visMax, vM, vMMma, vNMmpa, vNSigm, vP, vPhi, vS, vTheta, vTimes, vtl, vX1, vX2, vX3, vX4, vY1, vY2, vY3, vY4, & & x3, x4, xDatum, xMean, xSum, & & y1, y2, y3, y4, yDatum, yMean, ySum, & & wedge, zBAsth, zErr !--------------------------------------------------------------------- ! BEGINNING OF EXECUTABLE CODE ! *** KLUDGE ALERT ************************************************* ! Conversion of parameters (constants) to variables should logically ! have no effect, but in fact helped to suppress some spurious ! messages from the old IBM vS-FORTRAN compiler: mxNode = maxNod mxDOF = 2 * mxNode mxEl = maxEl mxFEl = maxFEl mxBn = maxBN mxStar = maxAtP mxGeo = maxGeo mxStr = maxStr mxSKS = maxSKS mxGSR = maxDat mxMOR = maxMOR mxAdj = maxAdj ! ****************************************************************** ! Mark most pLates as LACKING extensive attached slabs... DO i = 1, nPlate slab_q(i) = .FALSE. END DO ! ...except for these particular cases: slab_q( 8) = .TRUE. ! 8 = AU = Australia slab_q(14) = .TRUE. ! 14 = CL = Caroline slab_q(15) = .TRUE. ! 15 = CO = Cocos slab_q(21) = .TRUE. ! 21 = IN = India slab_q(22) = .TRUE. ! 22 = JF = Juan de Fuca slab_q(34) = .TRUE. ! 34 = NZ = Nazca slab_q(37) = .TRUE. ! 37 = PA = Pacific slab_q(39) = .TRUE. ! 39 = PS = Philippine Sea slab_q(40) = .TRUE. ! 40 = RI = Rivera slab_q(46) = .TRUE. ! 46 = SS = Solomon Sea ! N.B. This array is used when scoring upper-mantle seismic- ! anisotropy data, because in the Earth5 set of models, ! plates with slab_q(iPlate) = .TRUE. {those with ! extensive attached slabs} are assumed to have ! minimal and/or unresolvable basal shear tractions, ! and are excluded from the scoring process. ! That is, if total basal-strength traction for the ! Pacific plate is toward the Northwest, this might ! be primarily due to the attached slabs, and it does ! not preclude weak, distributed basal shear tractions ! that might point in a completely different direction! ! Write header on output file: WRITE (iUnitT, 1) 1 FORMAT ( / & &' =============================================================='/ & &' I PROGRAM -OrbScore2- I'/ & &' I A spherical-Earth, thin-plate program for scoring I'/ & &' I time-averaged (anelastic) deformation solutions I'/ & &' I computed by finite element program -Shells-. I'/ & &' I I'/ & &' I by I'/ & &' I Peter Bird I'/ & &' I Department of Earth, Planetary, and Space Sciences I'/ & &' I University of California I'/ & &' I Los Angeles, CA 90095-1567 I'/ & &' I Version of 14 February 2021 I'/ & &' ==============================================================') wedge = ABS(90.0D0 - ABS(dipMax)) * radians_per_degree ! converted to radians away from vertical slide = subDip * radians_per_degree ! converted to radians ! --------------------------------------------------------------------- ! Preview .FEG file to determine array sizes: WRITE (*, 101) iUnitG 101 FORMAT (/' Attempting to read finite element grid from unit', I3/) READ (iUnitG, * , IOSTAT = ios) IF (ios /= 0) THEN WRITE(*, "(' ERR','OR: File not found, or file is empty,' & & /' or file is too short.')") CALL Pause() STOP END IF READ (iUnitG, * , IOSTAT = ios) numNod IF (ios /= 0) THEN WRITE(*, "(' ERR','OR: File not found, or file is empty,' & & /' or file is too short.')") CALL Pause() STOP END IF IF (numNod > mxNode) THEN WRITE (*, "(' ERR','OR: Increase PARAMETER maxNode to at least ', I8, ' and recompile OrbScore2.')") numNod CALL Pause() STOP END IF DO 102 i = 1, numNod READ (iUnitG, * , IOSTAT = ios) IF (ios /= 0) THEN WRITE(*, "(' ERR','OR:File not found, or file is empty,' & & /' or file is too short.')") CALL Pause() STOP END IF 102 CONTINUE READ (iUnitG, * , IOSTAT = ios) numEl IF (ios /= 0) THEN WRITE(*, "(' ERR','OR: File not found, or file is empty,' & & /' or file is too short.')") CALL Pause() STOP END IF IF (numeL > mxEl) THEN WRITE (*, "(' ERR','OR: Increase PARAMETER maxEl to at least ', I8, ' and recompile OrbScore2.')") numEl CALL Pause() STOP END IF !Initialize survey to find LRn = MAX(continuum_LRi(1:mxEl), fault_LRi(1:MXFel) LRn = 0 ! until incremented below... DO 103 i = 1, numEl READ (iUnitG, "(A)", IOSTAT = ios) longer_line IF (ios /= 0) THEN WRITE(*, "(' ERR','OR:File not found, or file is empty,' & & /' or file is too short.')") CALL Pause() STOP END IF CALL Extract_LRi (longer_line, & ! input & LRi, shorter_line) ! output LRn = MAX(LRn, LRi) 103 CONTINUE nFl = 0 READ (iUnitG, * , IOSTAT = ios) n IF (ios == 0) nFl = n nFl = MAX(nFl, 0) IF (nFl > mxFEl) THEN WRITE (*, "(' ERR','OR: Increase PARAMETER maxFEl to at least ', I8, ' and recompile OrbScore2.')") nFl CALL Pause() STOP END IF DO 105 i = 1, nFl READ (iUnitG, "(A)", IOSTAT = ios) longer_line IF (ios /= 0) THEN WRITE(*, "(' ERR','OR:File not found, or file is empty,' & & /' or file is too short.')") CALL Pause() STOP END IF CALL Extract_LRi (longer_line, & ! input & LRi, shorter_line) ! output LRn = MAX(LRn, LRi) 105 CONTINUE REWIND (UNIT = iUnitG) ! to prepare for CALL GetNet, below... ! DIMENSIONs using newly-determined size variable LRn: ALLOCATE ( LR_is_defined(0:LRn) ) ALLOCATE ( LR_is_used(0:LRn) ) LR_is_defined = .FALSE. ! whole array, until information is read, below... LR_is_used = .FALSE. ! whole array, until information is read, below... ALLOCATE ( LR_set_fFric(0:LRn) ) ALLOCATE ( LR_set_cFric(0:LRn) ) ALLOCATE ( LR_set_Biot(0:LRn) ) ALLOCATE ( LR_set_Byerly(0:LRn) ) ALLOCATE ( LR_set_aCreep(1:2, 0:LRn) ) ALLOCATE ( LR_set_bCreep(1:2, 0:LRn) ) ALLOCATE ( LR_set_cCreep(1:2, 0:LRn) ) ALLOCATE ( LR_set_dCreep(1:2, 0:LRn) ) ALLOCATE ( LR_set_eCreep(0:LRn) ) !Just for ease in debugging, initialize all (currently) undefined array values as zero: LR_set_fFric = 0.0D0 LR_set_cFric = 0.0D0 LR_set_Biot = 0.0D0 LR_set_Byerly = 0.0D0 LR_set_aCreep = 0.0D0 LR_set_bCreep = 0.0D0 LR_set_cCreep = 0.0D0 LR_set_dCreep = 0.0D0 LR_set_eCreep = 0.0D0 ! --------------------------------------------------------------------- ! Input finite element grid and data values at node points: CALL GetNet (iUnitG, iUnitT, & ! input & mxDOF, mxEl, mxFEl, mxNode, & & brief, continuum_LRi, cooling_curvature, & ! output & density_anomaly, & & dQdTdA, elev, fault_LRi, fDip, & & nFakeN, nFl, nodeF, nodes, nRealN, & & numEl, numNod, n1000, offMax, offset, & & title1, tLNode, xNode, yNode, zMNode, & & checkE, checkF, checkN) ! work ! Read scalar parameters: CALL ReadPm (iUnitP, iUnitT, names, nPlate, offMax, & ! input & alphaT, conduc, & ! output & d_fFric, d_cFric, d_Biot, d_Byerly, d_aCreep, d_bCreep, d_cCreep, d_dCreep, d_eCreep, & & everyP, gMean , & & gradie, iConve, iPVRef, & & maxItr, OKDelV, OKToQt, oneKm, radio , & & radius, refStr, rhoAst, rhoBar, rhoH2O, & & tAdiab, tauMax, temLim, title3, & & trHMax, tSurf , vTimes, zBAsth) ! --------------------------------------------------------------------- !Remember the default ("d_") Lithospheric Rheology as LR0, or LR_set_XXXX(0): LR_set_fFric(0) = d_fFric LR_set_cFric(0) = d_cFric LR_set_Biot(0) = d_Biot LR_set_Byerly(0) = d_Byerly LR_set_aCreep(1:2, 0) = d_aCreep(1:2) LR_set_bCreep(1:2, 0) = d_bCreep(1:2) LR_set_cCreep(1:2, 0) = d_cCreep(1:2) LR_set_dCreep(1:2, 0) = d_dCreep(1:2) LR_set_eCreep(0) = d_eCreep LR_is_defined(0) = .TRUE. ! Obtain extra input file with Lithospheric Rheologies from the user: IF (LRn > 0) THEN WRITE (*, *) WRITE (*, "(' Lithospheric Rheology indeces from 0 to ', I8, ' are used in this .feg file.')") LRn WRITE (*, 113) iUnitLR 113 FORMAT (/' Attempting to read table of needed Lithospheric Rheologies from unit', I3/) READ (iUnitLR, * , IOSTAT = ios) ! READ (and discard) column-header line at top IF (ios /= 0) THEN WRITE(*, "(' ERR', 'OR: File not found, or file is empty,' / ' or file is too short.')") CALL Pause() STOP END IF collect_LRs: DO READ (iUnitLR, *, IOSTAT = ios) i IF (ios /= 0) EXIT collect_LRs ! at EOF, probably IF ((i < 1).OR.(i > LRn)) THEN WRITE (*, "(' ERR', 'OR: LR#', I8, ' is outside the legal range of (1:', I8, ').')") i, LRn WRITE (*, "(' To make it legal, some element in the .feg file must use this (or higher) LR#.')") CALL Pause() STOP END IF BACKSPACE(iUnitLR) READ (iUnitLR, *, IOSTAT = ios) i, LR_set_fFric(i), LR_set_cFric(i), LR_set_Biot(i), LR_set_Byerly(i), & & LR_set_aCreep(1:2, i), LR_set_bCreep(1:2, i), LR_set_cCreep(1:2, i), LR_set_dCreep(1:2, i), & & LR_set_eCreep(i) IF (ios == 0) THEN LR_is_defined(i) = .TRUE. ELSE WRITE (*, "(' ERR', 'OR while trying to read 13 REAL*8 values that make up LR#', I8)") i CALL Pause() STOP END IF END DO collect_LRs CLOSE (iUnitLR) !Now, "stress-test" the continuum elements to be sure that each has a defined rheology: DO j = 1, numEl i = continuum_LRi(j) IF (.NOT.LR_is_defined(i)) THEN WRITE (*, "(' ERR', 'OR: Continuum element ', I8,' uses LR#', I8, ' which has NOT been defined!')") j, i CALL Pause() STOP ELSE LR_is_used(i) = .TRUE. END IF END DO !Now, "stress-test" the fault elements to be sure that each has a defined rheology: IF (nFl > 0) THEN DO j = 1, nFl i = fault_LRi(j) IF (.NOT.LR_is_defined(i)) THEN WRITE (*, "(' ERR', 'OR: Fault element ', I8,' uses LR#', I8, ' which has NOT been defined!')") j, i CALL Pause() STOP ELSE LR_is_used(i) = .TRUE. END IF END DO END IF !Write a report to the log-file, to provide a record of the LRs used: WRITE (iUnitT, *) WRITE (iUnitT, "(' ===========================================================================================================================')") WRITE (iUnitT, "(' Table of alternative Lithospheric Rheologies defined and used:')") WRITE (iUnitT, "(' LR# fFric cFric Biot Byerly aCreep(1) aCreep(2) bCreep(1) bCreep(2) cCreep(1) cCreep(2) dCreep(1) dCreep(2) eCreep')") DO i = 0, LRn IF (LR_is_defined(i).AND.LR_is_used(i)) THEN WRITE (iUnitT, "(' ', I8, F6.3, F6.3, F6.3, F7.3, ES10.2, ES10.2, F10.0, F10.0, F10.4, F10.4, ES10.2, ES10.2, F10.5)") & & i, LR_set_fFric(i), LR_set_cFric(i), LR_set_Biot(i), LR_set_Byerly(i), & & LR_set_aCreep(1:2, i), LR_set_bCreep(1:2, i), LR_set_cCreep(1:2, i), LR_set_dCreep(1:2, i), & & LR_set_eCreep(i) END IF END DO WRITE (iUnitT, "('===========================================================================================================================')") WRITE (iUnitT, *) END IF ! LRn > 0 ! --------------------------------------------------------------------- ! Check topology, and compute geometric properties: log_strike_adjustments = .FALSE. skipBC = .TRUE. CALL Square (brief, fDip, iUnitT, & & log_strike_adjustments, & & mxBn, mxEl, mxFEl, mxNode, & & mxStar, nFl, nodeF, nodes, & & numEl, numNod, skipBC, radius, wedge, & ! INTENT(IN) & xNode, yNode, & ! INTENT(INOUT) & area, detJ, & & dxs, dys, dxsp, dysp, edgeFS, & & edgeTs, fLen, fPFlt, fpsfer, & & fArg, nCond, nodCon, sita, & ! INTENT(OUT) & checkN, list) ! WORKSPACE sphere = (nCond == 0) ! Input velocities at node points: WRITE(iUnitT, 10) iUnitV 10 FORMAT(/' ATTEMPTING TO READ VELOCITIES OF NODES ON UNIT', I3/) CALL OldVel (iUnitT, iUnitV, mxNode, numNod, & ! INTENT(IN) & haveNV, title1, title2, title3, v) ! INTENT(OUT) IF (.NOT.haveNV) THEN WRITE (iUnitT, 20) iUnitV 20 FORMAT (/' ERROR: No nodal velocities found on unit ', I2) CALL Pause() STOP END IF !---------------------------------------------------------------- ! COMPUTE RMS VALUE OF VELOCITIES OF NODES ! (Perhaps RMS value of surface integral would be more ! elegant, but this is much easier and almost the same.) ! First, it is necessary to remove any net rotation: numAdj = numNod IF (mxAdj < numAdj) THEN WRITE (iUnitT, 220) numNod 220 FORMAT (/' In order to remove net rotation from model velocities' & & /' (so as to accurately assess RMS velocity),' & & /' it is necessary to send a dummy set of (0, 0)' & & /' benchmark velocities to subprogram Adjust.' & & /' This in turn requires that parameter maxAdj' & & /' must be at least equal to the number of' & & /' nodes, which currently is: ', I7, '. Please fix!') CALL Pause() STOP END IF DO 250 i = 1, numNod nLat(i) = 90.0D0 - xNode(i) * oezOPi eLon(i) = yNode(i) * oezOPi predic(1, i) = v(1, i) predic(2, i) = v(2, i) data(1, i) = 0.0D0 data(2, i) = 0.0D0 250 CONTINUE IF (floats) THEN ALLOCATE ( weights(numNod) ) weights = 1.0D0 ! whole array WRITE (iUnitT, 251) 251 FORMAT (/' Calling Adjust with equal weight on each node:') CALL Adjust (data, eLon, iUnitT, mxAdj, numNod, nLat, & ! INTENT(IN) & radius, weights, & & predic, & ! INTENT(INOUT) & eLonP, nLatP, rate) ! INTENT(OUT) DEALLOCATE (weights) END IF sum = 0.0D0 DO 290 i = 1, numNod sum = sum + predic(1, i)**2 + predic(2, i)**2 290 CONTINUE rmsv = DSQRT(sum / numNod) t2 = 1000.0D0 * secPYr * rmsv IF (floats) THEN rate = -rate t = ABS(rate) * 1000.0D0 * secPYr * radius WRITE (iUnitT, 295) rate, t, eLonP, nLatP, rmsv, t2 295 FORMAT (/' After removal of net rotation of ', ES10.2, & & ' radians/second (', F5.1, ' mm/yr at equator)' & & /' about a pole at ', F7.2, ' degrees East, ', F6.2, & & ' degrees North,' & & /' the RMS nodal velocity for this model is ', & & ES10.2, & & ' (', F7.2, ' mm/yr).') ELSE WRITE (iUnitT, 296) rmsv, t2 296 FORMAT (/' The RMS nodal velocity for this model is ', & & ES10.2, & & ' (', F7.2, ' mm/yr).') END IF !------------------------------------------------------------------- ! ********* GEODESY SCORING SECTION ********** ! Read geodetic data for scoring: ! (Dislocation-in-elastic-halfspace corrections for ! temporary locking of brittle parts of faults within ! the model were added in Summer 2000 by Zhen Liu and Peter Bird.) ! Note that this does NOT include corrections for ! temporary locking of any subduction zones which might lie ! along the boundary of a model, UNLESS fault elements ! are used to represent those subduction shear zones.) CALL GetGEO (iUnitB, iUnitT, iUnitY, iUnitC, mxGeo, & & pltGeo, pltCos, pltInt, & ! INTENT(IN) & geoTag, geoPhi, geoThe, & ! INTENT(OUT) & geoVel, geoSig, geoAzi, & & gpsFMT, numGeo) ! Flag "explod" indicates that number of benchmarks ! is equal to number of nodes. This happens only when ! the "geodesy" (.gps) file was created by Explode3, and ! has "benchmark" locations for every node in the original ! .feg file, but slightly displaced away from nodes. ! In this case, processing is identical, but an extra ! output file is produced to facilitate plotting ! the smoothed (interseismic) velocities. explod = (numGeo == numNod) IF (numGeo <= 0) THEN geodes = 0.0D0 CLOSE (UNIT = iUnitY, DISP = "DELETE") GO TO 400 END IF WRITE (iUnitT, 301) numGeo, iUnitB 301 FORMAT (/' ', I6, ' Geodetic velocity data were read from unit ', I2) !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - ! Find weights for geodetic benchmarks, based on associated nearest-neighbor surface areas: IF (DIs_Initialization_Needed()) THEN subdivision = 4 WRITE (iUnitT, "(/' Initializing uniform global grid at subdivision ', & & I1,' for area-weighting of data...')") subdivision CALL DInitialize_Weighting (subdivision) END IF WRITE (iUnitT, *) WRITE (iUnitT, "(' Computing area-weights for geodetic benchmarks...')") ALLOCATE ( weights(numGeo) ) CALL DPerform_Weighting (number_of_data = numGeo, & & theta_radians = geoThe, & & phi_radians = geoPhi, & ! INTENT(IN) & weights = weights) ! INTENT(OUT) !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - ! Find predicted velocities at benchmarks ! (A) Steady-state components: IF (explod) THEN WRITE (iUnitT, 302) 302 FORMAT(/' Calling FindIt to determine element containing' & & /' each displaced node; be patient!....') END IF i = 0 303 i = i + 1 theta = geoThe(i) phi = geoPhi(i) CALL FindIt (mxEl, mxNode, nodes, numEl, xNode, yNode, & ! INTENT(IN) & theta, phi, & & iEle, s1, s2, s3) ! INTENT(OUT) IF (iEle == 0) THEN ! Bemchmark is outside of grid; delete it. WRITE (iUnitT, 305) geoTag(i) 305 FORMAT (' Benchmark ', A20, ' was outside the grid: deleted.') DO 310 j = i, numGeo - 1 geoTag(j) = geoTag(j + 1) geoPhi(j) = geoPhi(j + 1) geoThe(j) = geoThe(j + 1) geoVel(j) = geoVel(j + 1) geoSig(j) = geoSig(j + 1) geoAzi(j) = geoAzi(j + 1) 310 CONTINUE numGeo = numGeo - 1 i = i - 1 ELSE theLat = 90.0D0 - theta * oezOPi theLon = phi * oezOPi IF (theLon > +180.0D0) theLon = theLon - 360.0D0 IF (theLon > +180.0D0) theLon = theLon - 360.0D0 IF (theLon < -180.0D0) theLon = theLon + 360.0D0 IF (theLon < -180.0D0) theLon = theLon + 360.0D0 IF (.NOT.explod) THEN showis = .FALSE. IF (showis) THEN WRITE(iUnitT, 315)geoTag(i), theLon, theLat, iEle, s1, s2, s3 315 FORMAT(' ', A20,' (', F7.2, 'E, ', F6.2, 'N) = #', I5, ':(', F5.3, ',', F5.3, ',', F5.3,')') END IF END IF ! Interpolate v at location (s1, s2, s3) in element #iEle: CALL Projec (iEle, & & mxEl, mxNode, nodes, & & s1, s2, s3, & & xNode, yNode, v, & ! INTENT(IN) & vTheta, vPhi) ! INTENT(OUT) predic(1, i) = vTheta predic(2, i) = vPhi longterm_at_GPS(1, i) = vTheta ! This duplicate copy of the long-term velocities longterm_at_GPS(2, i) = vPhi ! at GPS benchmarks (predicted by Shells) will NOT ! be modified by CALL Adjust (or in any other way), ! so that it will remain available IF (pltInt). data(1, i) = geoVel(i) * COS(geoAzi(i)) data(2, i) = geoVel(i) * SIN(geoAzi(i)) nLat(i) = 90.0D0 - geoThe(i) * oezOPi eLon(i) = geoPhi(i) * oezOPi END IF ! End of variable-length loop: IF (i < numGeo) GO TO 303 IF (explod.AND.(numGeo < numNod)) THEN ! Condition (numGeo == numNod) that made "explod" .TRUE. is no longer valid: explod = .FALSE. WRITE (iUnitT, 319) 319 FORMAT(/' This should not happen OFTEN when benchmark' & & /' locations are really node locations, slightly' & & /' displaced by Explode3.' & & /' Probably the bad benchmark(s) is/are only' & & /' slightly outside the domain of the F-E model.' & & /" Use OrbWin (or OrbWeaver) to view the Explode3'd" & & /' grid and see why some nodes/benchmarks moved' & & /' out of bounds. Then adjust the corresponding' & & /' benchmark coordinates in the false-benchmark' & & /' data file, and re-run OrbScore2.') CALL Pause () STOP END IF !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - ! Find predicted velocities at benchmarks; ! (B) Coseismic displacement rate due to brittle ! parts of faults within the model domain. ! Find thicknesses of crust and mantle-lithosphere: CALL Interp (zMNode, mxEl, mxNode, nodes, numEl, & ! INTENT(IN) & zMoho) ! INTENT(OUT) CALL Interp (tLNode, mxEl, mxNode, nodes, numEl, & ! INTENT(IN) & tLInt) ! INTENT(OUT) ! Compute tactical values of limits on viscosity, and weights for ! imposition of constraints in linear systems: CALL Limits (area, detJ, iUnitT, mxEl, numEl, & & OKDelV, radius, refStr, sphere, tLInt, & & trHMax, zMoho, & ! INTENT(IN) & constr, etaMax, fMuMax, visMax) ! INTENT(OUT) ! Locate the brittle/ductile transition in each fault: ! Begin with very crude initialization: DO i = 1, nFl zTranF(1, i) = MAX(zMNode(nodeF(1, i)) / 2.0D0, oneKm) zTranF(2, i) = MAX(tLNode(nodeF(1, i)) / 2.0D0, oneKm) END DO ! Then search for actual brittle/ductile transition, iteratively: DO iMohr = 1, 3 CALL Mohr (alphaT, conduc, constr, & ! input & continuum_LRi, & & dQdTdA, elev, & & fault_LRi, fDip, fMuMax, & & fPFlt, fArg, gMean, & & LRn, LR_set_fFric, LR_set_cFric, LR_set_Biot, LR_set_Byerly, & & LR_set_aCreep, LR_set_bCreep, LR_set_cCreep, LR_set_dCreep, LR_set_eCreep, & & mxEl, mxFEl, mxNode, nFl, nodeF, & & offMax, offset, & & oneKm, radio, rhoH2O, rhoBar, & & slide, tauMax, & & tLNode, tSurf, v, wedge, & & zMNode, & & zTranF, & ! modify & fC, fIMuDZ, fPeakS, fSlips, fTStar) ! output END DO ! NEW FEATURE (2018.11.26): ! In all subduction zones (oceanic or continental), ! override the brittle/ductile transition(s) ! determined (above) by Mohr, and replace them with a single 26-km-thick ! frictional layer, following Oleskevich et al. [1999, J. Geophys. Res.] ! and Bird & Kagan [2004, Bull. Seismol. Soc. Am.]: DO i = 1, nFl dip = 0.5D0 * (fDip(1, i) + fDip(2, i)) ! averaging the dips at the 2 ends of fault trace. IF ((dip <= slide).OR.((Pi - dip) <= slide)) THEN ! N.B. slide is limiting (maximum) value for SUBs, from PARAMETER subDip, after conversion to radians. !This fault is a SUB: zTranF(1, i) = 26.0D3 ! 26 km = average depth range of seismogenic patches in SUBs, from Oleskevich et al. [1999, JGR]. zTranF(2, i) = 0.0D3 ! (Not adding any other frictional layer at the top of the mantle, which is probably serpentinized.) END IF END DO ! Compute benchmark velocities due to elastic dislocations ! from earthquakes in the brittle part of each fault, at the slip-rate ! predicted by the -Shells- finite-element model: IF (explod) THEN WRITE (iUnitT, 340) 340 FORMAT(/' Calling Coseis to evaluate long-term-average coseismic' & & /' velocity of each node; be patient!....') ELSE WRITE (iUnitT, 341) 341 FORMAT(/' Calling Coseis to evaluate long-term-average coseismic' & & /' velocity of each benchmark....') END IF CALL DCoseis (fArg, fDip, geoThe, geoPhi, mxNode, & & mxFEl, numGeo, nFl, nodeF, radius, & & v, wedge, xNode, yNode, zMNode, zTranF, & ! INTENT(IN) & geoUTh, geoUPh) ! INTENT(OUT) ! For a typical year with no earthquakes, ! combine the predicted benchmark velocities ! by subtracting the coseismic part: DO 350 i = 1, numGeo predic(1, i) = predic(1, i) - geoUTh(i) predic(2, i) = predic(2, i) - geoUPh(i) 350 CONTINUE IF (explod) THEN WRITE (iUnitT, 351) 351 FORMAT(//' Creating file with smoothed (interseismic)' & & /' velocities in v_____.out format, to be used' & & /' as input to FiniteMap for visualizing them.') WRITE(iUnitT, 352) iUnitX 352 FORMAT(/' Attempting to create output v____.out file for' & & /' smoothed (interseismic) velocities at benchmarks' & & /' (which were slightly displaced from nodes by' & & /' Explode3) using unit ',I3/) OPEN (UNIT = iUnitX, FILE = '') line = "Explode3 of : "//TRIM(title1(1:66)) WRITE(iUnitX, "(A)") TRIM(line) line = "Interseismic: "//TRIM(title2(1:66)) WRITE(iUnitX, "(A)") TRIM(line) line = "velocities : "//TRIM(title3(1:66)) WRITE(iUnitX, "(A)") TRIM(line) WRITE(iUnitX, "(4ES18.9)")((predic(i, j), i = 1, 2), j = 1, numGeo) ! Note: numGeo == numNod, or we wouldn't be here. CLOSE(iUnitX) WRITE (iUnitT, 353) 353 FORMAT(/' Smoothed velocities have been written.' & & /' This run of OrbScore2 will now end, as' & & /' it was not intended for actual scoring,' & & /' but just to evaluate interseismic' & & /' velocities at ""benchmarks"" which are' & & /' really slightly-displaced nodes.') CALL Pause () STOP END IF !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - ! Adjust predictions by rigid rotation to best-fit geodesy: WRITE (iUnitT, 362) 362 FORMAT (/' Calling Adjust with area-weight on each benchmark' & & /' to minimize global area-integral of velocity-error-squared:') IF (floats) THEN CALL Adjust (data, eLon, iUnitT, mxAdj, numGeo, nLat, & & radius, weights, & ! INTENT(IN) & predic, & ! INTENT(INOUT) & eLonP, nLatP, rate) ! INTENT(OUT) t = ABS(rate) * 1000.0D0 * secPYr * radius WRITE (iUnitT, 370) rate, t, eLonP, nLatP 370 FORMAT (/' Model predictions adjusted by rotation of ', ES10.2, & & ' radians/second (', F7.2,' mm/yr at equator)' & & /' about a pole at ', F7.2, ' degrees East, ', F6.2, & & ' degrees North,' & & /' before comparison to geodetic data.') END IF WRITE (iUnitT, 380) 380 FORMAT (// & &' Geodetic Benchmark Velocities Versus Model Predictions:'// & &' Benchmark Velocity (mm/yr) Azimuth Model.v (mm/yr)' & &,' Model.Az Error (mm/yr) (Sigmas) area-weight'/ & &' --------- -------- ----- ------- ------- -----', & &' -------- ----- ----- ------ ----------') 381 FORMAT (' ', A20, ES9.2, ' (', F5.1, ') ', F7.1, ES10.2, & &' (', F5.1, ') ', F8.1, ES10.2, ' (', F5.1, ') (', F5.1, ') ', F10.5) sum0s = 0.0D0 sum1 = 0.0D0 sum1s = 0.0D0 sum2 = 0.0D0 sum2s = 0.0D0 sumN = 0.0D0 sumNs = 0.0D0 DO 390 i = 1, numGeo gVMma = geoVel(i) * secPYr * 1000.0D0 azim = 180.0D0 - geoAzi(i) * oezOPi IF (azim < 0.0D0) azim = azim + 360.0D0 IF (azim > 360.0D0) azim = azim - 360.0D0 vM = DSQRT(predic(1, i)**2 + predic(2, i)**2) vMMma = vM * secPYr * 1000.0D0 pAzim = 180.0D0 - ATan2F(predic(2, i), predic(1, i)) * oezOPi IF (pAzim < 0.0D0) pAzim = pAzim + 360.0D0 IF (pAzim > 360.0D0) pAzim = pAzim - 360.0D0 bad = DSQRT((predic(1, i) - data(1, i))**2 + & & (predic(2, i) - data(2, i))**2) badMma = bad * secPYr * 1000.0D0 badSig = bad / geoSig(i) WRITE (iUnitT, 381) geoTag(i), geoVel(i), gVMma, azim, & & vM, vMMma, pAzim, bad, badMma, badSig, & & weights(i) IF (pltGeo) THEN eLonDe = geoPhi(i) / piO180 nLatDe = 90.0D0 - (geoThe(i) / piO180) nLatDe = MIN(90.0D0, MAX(-90.0D0, nLatDe)) vEMmpa = +(predic(2, i) - data(2, i)) * 1000.0D0 * secPYr ! comparing (reference-frame-adjusted Shells predictions vNMmpa = -(predic(1, i) - data(1, i)) * 1000.0D0 * secPYr ! of interseismic velocities) at benchmarks to (data) vESigm = 0.1D0 vNSigm = 0.1D0 correl = 0.0D0 frame = "[same]" tag = geoTag(i) WRITE (iUnitY, gpsFMT) eLonDe, nLatDe, vEMmpa, vNMmpa, & & vESigm, vNSigm, correl, & & TRIM(frame), TRIM(tag) END IF ! pltGeo IF (pltCos) THEN ! write out the Coseismic part of the Shells-model predicted geodetic velocities: eLonDe = geoPhi(i) / piO180 nLatDe = 90.0D0 - (geoThe(i) / piO180) nLatDe = MIN(90.0D0, MAX(-90.0D0, nLatDe)) vEMmpa = +geoUPh(i) * 1000.0D0 * secPYr vNMmpa = -geoUTh(i) * 1000.0D0 * secPYr vESigm = 0.1D0 vNSigm = 0.1D0 correl = 0.0D0 frame = "[same]" tag = geoTag(i) WRITE (iUnitC, gpsFMT) eLonDe, nLatDe, vEMmpa, vNMmpa, & & vESigm, vNSigm, correl, & & TRIM(frame), TRIM(tag) END IF ! pltCos IF (pltInt) THEN ! write out the (NON-frame-adjusted) Interseismic part of the Shells-model predicted geodetic velocities: eLonDe = geoPhi(i) / piO180 nLatDe = 90.0D0 - (geoThe(i) / piO180) nLatDe = MIN(90.0D0, MAX(-90.0D0, nLatDe)) vEMmpa = +(longterm_at_GPS(2, i) - geoUPh(i)) * 1000.0D0 * secPYr ! long-term Shells velocity minus mean co-seismic velocity, vNMmpa = -(longterm_at_GPS(1, i) - geoUTh(i)) * 1000.0D0 * secPYr ! both evaluated at geodetic benchmarks vESigm = 0.1D0 vNSigm = 0.1D0 correl = 0.0D0 frame = "[same]" tag = geoTag(i) WRITE (iUnitI, gpsFMT) eLonDe, nLatDe, vEMmpa, vNMmpa, & & vESigm, vNSigm, correl, & & TRIM(frame), TRIM(tag) END IF ! pltInt IF (badSig > 2.0D0) sum0s = sum0s + weights(i) sum1 = sum1 + weights(i) * bad sum1s = sum1s + weights(i) * badSig sum2 = sum2 + weights(i) * bad**2 sum2s = sum2s + weights(i) * badSig**2 sumN = MAX(sumN, bad) sumNs = MAX(sumNs, badSig) 390 CONTINUE sum0s = sum0s / numGeo sum1 = sum1 / numGeo sum1mm = sum1 * 1000.0D0 * secPYr sum1s = sum1s / numGeo sum2 = SQRT(sum2 / numGeo) sum2mm = sum2 * secPYr * 1000.0D0 sum2s = SQRT(sum2s / numGeo) sumNmm = sumN * 1000.0D0 * secPYr WRITE (iUnitT, 397) sum0s, sum1, sum1mm, sum1s, & & sum2, sum2mm 397 FORMAT(/ & & /' Summary of geodetic prediction errors:' & & /' Fraction of predictions wrong by over 2-sigma: ', F5.3 & & /' Mean velocity error: ', ES10.2, ' (', F7.2, ' mm/yr)' & & /' Mean number of sigmas in error: ', F5.2, & & /' RMS velocity error: ', ES10.2, ' (', F7.2, ' mm/yr)' & &) IF (floats) THEN WRITE (iUnitT, 398) 398 FORMAT(' *** Note: Line above was optimized by Adjust. ***') END IF WRITE (iUnitT, 399) sum2s, sumN, sumNmm, sumNs 399 FORMAT(' RMS number of sigmas in error: ', F5.2 & & /' Worst velocity error :', ES10.2, ' (', F7.2, ' mm/yr)' & & /' Largest error, in sigmas: ', F6.2) ! Select one indicator of error (RMS error, after frame-change): geodes = sum2mm IF (pltGeo) THEN CLOSE(iUnitY) WRITE(iUnitT, "(/' Wrote ERRORS.gps, for plotting geodetic errors with -FiniteMap-.')") END IF IF (pltCos) THEN CLOSE(iUnitC) WRITE(iUnitT, "(/' Wrote COSEISMIC.gps, for plotting model coseismic velocities with -FiniteMap-.')") END IF IF (pltInt) THEN CLOSE(iUnitI) WRITE(iUnitT, "(/' Wrote INTERSEISMIC.gps, for plotting model interseismic velocities with -FiniteMap-.')") END IF DEALLOCATE (weights) ! ********* END SPACE-GEODESY SCORING SECTION ********** !------------------------------------------------------------------- ! ****** STRESS-DIRECTION SCORING SECTION ******** ! Read most-compressive horizontal principal stress azimuths for scoring. 400 CALL GetSTR (iUnitS, iUnitT, mxStr, & ! INTENT(IN) & strTag, strThe, strPhi, & & strArg, strQua, strReg, numStr) ! INTENT(OUT) IF (numStr <= 0) THEN stress = 0.0D0 regime = 0.0D0 GO TO 1301 END IF WRITE (iUnitT, 401) numStr, iUnitS 401 FORMAT (/' ',I6,' Stress direction data were read from unit ', I2) !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - IF (DIs_Initialization_Needed()) THEN subdivision = 4 WRITE (*, "(/' Initializing uniform global grid at subdivision ', & & I1, ' for area-weighting of data...')") subdivision CALL DInitialize_Weighting (subdivision) END IF WRITE (*, "(' Computing area-weights for stress direction/regime data...')") ALLOCATE ( weights(numStr) ) weights = 1.0D0 ! whole array ! CALL DPerform_Weighting (number_of_data = numStr, & ! & theta_radians = strThe, & ! & phi_radians = strPhi, & ! INTENT(IN) ! & weights = weights) ! INTENT(OUT) ! N. B. Area-weighting turned off here when we switched to use of ! robust_interpolated_stress_for_OrbScore2, 2007.03.12 in the Earth5 project. !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ! Compute the (theta, phi) and then Cartesian (x, y, z) of integration points: DO 450 i = 1, numEl DO 420 k = 1, 3 CartVs(1, k) = DCOS(yNode(nodes(k, i))) * & & DSIN(xNode(nodes(k, i))) CartVs(2, k) = DSIN(yNode(nodes(k, i))) * & & DSIN(xNode(nodes(k, i))) CartVs(3, k) = DCOS(xNode(nodes(k, i))) 420 CONTINUE DO 440 m = 1, 7 DO 430 j = 1, 3 tempV(j) = 0.0D0 DO 425 k = 1, 3 tempV(j) = tempV(j) + CartVs(j, k) * points(k, m) 425 CONTINUE 430 CONTINUE CALL Unit(tempV) ! INTENT(INOUT) IF (ABS(tempV(3)) <= 0.5D0) THEN xIP(m, i) = DACOS(tempV(3)) ELSE equPar = DSQRT(tempV(1)**2 + tempV(2)**2) xIP(m, i) = ATan2F(equPar, tempV(3)) END IF yIP(m, i) = ATan2F(tempV(2), tempV(1)) CartR(1, m, i) = DSIN(xIP(m, i)) * DCOS(yIP(m, i)) CartR(2, m, i) = DSIN(xIP(m, i)) * DSIN(yIP(m, i)) CartR(3, m, i) = DCOS(xIP(m, i)) 440 CONTINUE 450 CONTINUE ! Strain-rates at integration points: CALL EDot (dxs, dys, & & fpsfer, mxEl, & & mxNode, nodes, numEl, radius, sita, v, & ! INTENT(IN) & eRate) ! INTENT(OUT) WRITE (iUnitT, 460) 460 FORMAT (// & &' Horizontal most-compressive stress directions versus ', & & 'model predictions:'/ & &/' Datum E.Long. N.Latt. Element Point', & & ' Quality Azimuth Model.Az Error', & & ' Regime E1 E2 E3 Err Match? area-weight' & &/' ----- ------- ------- ------- ----- ------- ', & & '------- -------- -----', & & ' ------ -------- -------- -------- -------- ------ -----------') 461 FORMAT (' ', A5, 1X, F7.2, 1X, F7.2, 1X, I7, 1X, I5, 1X, F7.0, 1X, & & F7.0, 1X, F8.0, 1X, F5.0, & & 5X, A2, 4ES9.1, 1X, A4, 2X, F10.5) !: Accumulator of bad stress regimes: nBadRe = 0 ! simple integer count of bad regimes aBadRe = 0.0D0 ! area-weighted error count aRegime_denominator = 0.0D0 !: Accumulators for general data (all qualities): sum1 = 0.0D0 sum1d = 0.0D0 sum2 = 0.0D0 sum2d = 0.0D0 !: Accumulators for A,B-quality data only: n_AB_data = 0 aSum1 = 0.0D0 aSum1D = 0.0D0 aSum2 = 0.0D0 aSum2D = 0.0D0 ! NOTE: Integration point must be within one-half element- ! width of datum, or datum is not used. ! Here, the typical element side is deterMined by ! averaging together all element sides, ! and expressed in radians: sumSid = 0.0D0 DO 471 i = 1, numEl DO 470 k = 1, 3 kp1 = 1 + MOD(k, 3) n1 = nodes(k, i) n2 = nodes(kp1, i) x1 = xNode(n1) x2 = xNode(n2) y1 = yNode(n1) y2 = yNode(n2) rA(1) = DSIN(x1) * DCOS(y1) rA(2) = DSIN(x1) * DSIN(y1) rA(3) = DCOS(x1) rB(1) = DSIN(x2) * DCOS(y2) rB(2) = DSIN(x2) * DSIN(y2) rB(3) = DCOS(x2) dot = rA(1) * rB(1) + rA(2) * rB(2) + rA(3) * rB(3) CALL DCross (rA, rB, & ! INTENT(IN) & cUvec) ! INTENT(OUT) croSiz = DSQRT(cUvec(1)**2 + cUvec(2)**2 + cUvec(3)**2) side = ATan2F(croSiz, dot) sumSid = sumSid + side 470 CONTINUE 471 CONTINUE toler = sumSid / (6.0D0 * numEl) numDen = numStr DO 490 i = 1, numStr ! Next stress direction datum: theLon = oezOPi * strPhi(i) IF (theLon > 180.0D0) theLon = theLon - 360.0D0 IF (theLon < -180.0D0) theLon = theLon + 360.0D0 theLat = 90.0D0 - oezOPi * strThe(i) rS(1) = DSIN(strThe(i)) * DCOS(strPhi(i)) rS(2) = DSIN(strThe(i)) * DSIN(strPhi(i)) rS(3) = DCOS(strThe(i)) ! Find closest integration point in the grid: r2min = 9.99D29 DO 478 m = 1, 7 DO 476 j = 1, numEl r2 = (rS(1) - CartR(1, m, j))**2 + & & (rS(2) - CartR(2, m, j))**2 + & & (rS(3) - CartR(3, m, j))**2 IF (r2 < r2min) THEN iEle = j mIP = m r2min = r2 END IF 476 CONTINUE 478 CONTINUE rMin = DSQRT(r2min) IF (rMin > toler) THEN WRITE (iUnitT, 480) strtag(i), theLon, theLat, toler 480 FORMAT (' ',A5,1X,F7.2,1X,F7.2, & & ' Stress datum is more than ', & & F6.4, ' radians from nearest' & & ,' integration point: Ignored.') numDen = numDen - 1 ELSE ! Get sigma_1H direction from E_1 direction: e11h = eRate(1, mIP, iEle) e22h = eRate(2, mIP, iEle) e12h = eRate(3, mIP, iEle) CALL Prince (e11h, e22h, e12h, & ! INTENT(IN) & e1h, e2h, u1x, u1y, u2x, u2y) ! INTENT(OUT) err = -(e11h + e22h) ! (or, = -(e1h + e2h) ) e1 = MIN(e1h, e2h, err) e3 = MAX(e1h, e2h, err) IF ((err > e1).AND.(err < e3)) THEN e2 = err ELSE IF ((e1h > e1).AND.(e1h < e3)) THEN e2 = e1h ELSE e2 = e2h END IF pAzim = 180.0D0 - oezOPi * ATan2F(u1y, u1x) IF (pAzim < 0.0D0) pAzim = pAzim + 180.0D0 IF (pAzim < 0.0D0) pAzim = pAzim + 180.0D0 IF (pAzim >= 180.0D0) pAzim = pAzim - 180.0D0 IF (pAzim >= 180.0D0) pAzim = pAzim - 180.0D0 azim = 180.0D0 - oezOPi * strArg(i) IF (azim < 0.0D0) azim = azim + 180.0D0 IF (azim < 0.0D0) azim = azim + 180.0D0 IF (azim >= 180.0D0) azim = azim - 180.0D0 IF (azim >= 180.0D0) azim = azim - 180.0D0 bad = ABS(pAzim - azim) IF (bad > 90.0D0) bad = ABS(bad - 180.0D0) IF (bad > 90.0D0) bad = ABS(bad - 180.0D0) IF (bad > 90.0D0) bad = ABS(bad - 180.0D0) IF (bad > 90.0D0) bad = ABS(bad - 180.0D0) sum1 = sum1 + bad * strQua(i) *weights(i) sum1d = sum1d + strQua(i) * weights(i) sum2 = sum2 + strQua(i) * bad**2 * weights(i) sum2d = sum2d + strQua(i) * weights(i) ! Maintain separate totals for A,B-quality only: IF (strQua(i) > 2.9D0) THEN n_AB_data = n_AB_data + 1 aSum1 = aSum1 + bad * weights(i) aSum1D = aSum1D + weights(i) aSum2 = aSum2 + bad**2 * weights(i) aSum2D = aSum2D + weights(i) END IF ! Check whether stress regime is correct: reg = strreg(i) outcom = 'OK? ' IF (reg == 'NF') THEN ! Normal faulting; ! e1 should be vertical IF (e1 /= err) THEN nBadRe = nBadRe + 1 aBadRe = aBadRe + weights(i) outcom = 'BAD ' ELSE outcom = 'GOOD' END IF ELSE IF (reg == 'SS') THEN ! Strike-slip faulting; ! e2 should be vertical IF (e2 /= err) THEN nBadRe = nBadRe + 1 aBadRe = aBadRe + weights(i) outcom = 'BAD ' ELSE outcom = 'GOOD' END IF ELSE IF (reg == 'TF') THEN ! Thrust-faulting; ! e3 should be vertical IF (e3 /= err) THEN nBadRe = nBadRe + 1 aBadRe = aBadRe + weights(i) outcom = 'BAD ' ELSE outcom = 'GOOD' END IF ELSE IF ((reg == 'ST').OR.(reg == 'TS')) THEN ! Transpression; ! err should be positive, and ! the two horizontal principal stresses ! should have opposite signs. IF ((err <= 0.0D0).OR.((e1h * e2h) >= 0.0D0)) THEN nBadRe = nBadRe + 1 aBadRe = aBadRe + weights(i) outcom = 'BAD ' ELSE outcom = 'GOOD' END IF ELSE IF ((reg == 'NS').OR.(reg == 'SN')) THEN ! Transtension; ! err should be negative, and ! the two horizontal principal stresses ! should have opposite signs. IF ((err >= 0.0D0).OR.((e1h * e2h) >= 0.0D0)) THEN nBadRe = nBadRe + 1 aBadRe = aBadRe + weights(i) outcom = 'BAD ' ELSE outcom = 'GOOD' END IF END IF aRegime_denominator = aRegime_denominator + weights(i) WRITE (iUnitT, 461) strtag(i), theLon, theLat, iEle, mIP, & & strQua(i), azim, pAzim, bad, & & strreg(i), e1, e2, e3, ERR, outcom, weights(i) END IF 490 CONTINUE IF (sum1d > 0.0D0) THEN sum1 = sum1 / sum1d ELSE sum1 = 0.0D0 END IF IF (sum2d > 0.0D0) THEN sum2 = DSQRT(sum2 / sum2d) ELSE sum2 = 0.0D0 END IF WRITE (iUnitT, 496) numDen, sum1, sum2 496 FORMAT (/' Summary of Stress-Direction Errors:' & & /' Number of stress directions used for scoring: ', I6, '.' & & /' Area- & quality-weighted MEAN error: ', F5.2, ' degrees.' & & /' Area- & quality-weighted RMS error: ', F5.2, ' degrees.') IF ((aSum1D > 0.0D0).AND.(aSum2D > 0.0D0)) THEN aSum1 = aSum1 / aSum1D aSum2 = DSQRT(aSum2 / aSum2D) WRITE (iUnitT, 497) n_AB_data, aSum1, aSum2 497 FORMAT (/' OR, using only the ',I6,' A,B-quality data:' & & /' Area-weighted MEAN error: ', F5.2, ' degrees.' & & /' Area-weighted RMS error: ', F5.2, ' degrees.') END IF IF (numDen > 0.0D0) THEN perBad = (100.0D0 * nBadRe) / (1.0D0 * numDen) ELSE perBad = 0.0D0 END IF WRITE (iUnitT, 498) nBadRe, numDen, perBad 498 FORMAT (' Regime:', I6, ' bad stress regimes out of ', I6, & & ' data = ', F6.2, ' % bad.') IF (aRegime_denominator > 0.0D0) THEN perBad = (100.0D0 * aBadRe) / aRegime_denominator ELSE perBad = 0.0D0 END IF WRITE (iUnitT, 499) aBadRe, aRegime_denominator, perBad 499 FORMAT (' Regime:', F10.3, ' area-weighted bad stress regimes out of ', F10.3, & & ' area-weighted tests = ', F6.2, ' % bad.') ! Arbitrarily choose measure of error (mean error, in degrees): stress = sum1 DEALLOCATE ( weights ) ! ********* END STRESS SCORING SECTION ********** !------------------------------------------------------------------- ! ********* SLIP-RATE SCORING SECTION *********** ! Read geological slip rates on faults which are represented by ! fault elements in the .FEG file. ! Score both components: strike-slip, and relative-vertical (throw) data ! where available. ! This part of the code was modified from the -PLATES- scoring program, ! -Score_AK.for. ! The data format is given in: Slip_rate_format.txt. ! Read test dataset of geologic slip-rate limits for faults: 1301 CALL GetSLF(iUnitD, iUnitT, mxGSR, & ! INTENT(IN) & fName, rLat, delVZ, fltLon, fltLat, nFData) ! INTENT(OUT) IF (nFData <= 0) THEN slpErr = 0.0D0 GO TO 501 END IF ! Compute slip rates (in mm/year) at specified points ! and assign scores and probabilities: IF (table_of_offset_rates) THEN OPEN (UNIT = 51, FILE = "offset_rate_misfits.txt") END IF ! table_of_offset_rates (optional log file when running interactively) WRITE (iUnitT, 1310) 1310 FORMAT (/' Geologic sliprate predictions and their errors (in mm/year):') WRITE (iUnitT, 1320) 1320 FORMAT (' ' & & ' Minimum Maximum Model' & &,' Dextral Minimum Maximum Model Throw ' & &/' Fault name (first 100 bytes) Dextral Dextral Dextral' & &,' Error Chance Throw Throw Throw Error Chance' & &/' ---------------------------------------------------------------------------------------------------- ------- ------- -------' & &,' ------- ------ ------- ------- ------- ------- ------') IF (table_of_offset_rates) THEN WRITE (51, 1310) WRITE (51, 1320) END IF ! table_of_offset_rates (optional log file when running interactively) shrink = 1.0D0 nGSDat = 0 nB = 0 sumB = 0.0D0 sumB2 = 0.0D0 bigB = 0.0D0 DO 1390 j = 1, nFData ! Find fault element closest to (fltLon, fltLat): ! xDatum is Theta, S from North pole, in radians: xDatum = (90.0D0 - fltLat(j)) * piO180 ! yDatum is Phi, E from Greenwich meridan, in radians: yDatum = fltLon(j) * piO180 IF (nFl == 0) THEN WRITE(iUnitT, 1331) j 1331 FORMAT(' ERR0R: NO FAULT ELEMENTS IN THIS .FEG FILE' & & ' TO MATCH FAULT SLIP-RATE DATUM #', I8, '.') CALL Pause() STOP END IF i = 0 s = 0.0D0 r2min = 9.99D+29 DO 1340 k = 1, nFl n1 = nodeF(1, k) n2 = nodeF(2, k) n3 = nodeF(3, k) n4 = nodeF(4, k) x1 = xNode(n1) x2 = xNode(n2) x3 = xNode(n3) x4 = xNode(n4) y1 = yNode(n1) y2 = yNode(n2) y3 = yNode(n3) y4 = yNode(n4) factor = SIN(0.5D0 * (x1 + x2)) r2flt = (x2 - x1)**2 + (factor * (y1 - y2))**2 r2n1 = (xDatum - x1)**2 + (factor * (yDatum - y1))**2 r2n2 = (xDatum - x2)**2 + (factor * (yDatum - y2))**2 IF ((r2n1 <= r2flt).OR.(r2n2 <= r2flt)) THEN ! Determine internal coordinates for this fault: IF ((r2n1 <= r2flt).AND.(r2n2 <= r2flt)) THEN st = r2n1 / (r2n1 + r2n2) smallc = SQRT(r2flt) smallb = SQRT(r2n1) smalla = SQRT(r2n2) ! a = angle at the node-1 end: cosa = (r2n1 + r2flt - r2n2) / (2.0 * smallb * smallc) a = ACOS(cosa) ! b = angle at the node-2 end: cosb = (r2n2 + r2flt - r2n1) / (2.0 * smalla * smallc) b = ACOS(cosb) ! r = offset from the fault line: r = smallc * (SIN(a) * SIN(b) / SIN(a + b)) r2 = r * r ELSE IF (r2n1 <= r2flt) THEN st = 0.0 r2 = r2n1 ELSE IF (r2n2 <= r2flt) THEN st = 1.0 r2 = r2n2 END IF ! Replace current best solution? IF (r2 <= r2min) THEN i = k s = st r2min = r2 END IF END IF 1340 CONTINUE IF (i == 0) THEN WRITE(iUnitT, 1341) j, TRIM(fName(j)), fltLon(j), fltLat(j) 1341 FORMAT(' ERR0R: NO FAULT ELEMENT NEAR TO DATUM ', I4 & & /' ', A & & /' AT (', F10.4, ',', F10.4,')') CALL Pause() STOP END IF ! Assuming that i = fault-element was found, continue: n1 = nodeF(1, i) n2 = nodeF(2, i) n3 = nodeF(3, i) n4 = nodeF(4, i) x1 = xNode(n1) x2 = xNode(n2) x3 = xNode(n3) x4 = xNode(n4) y1 = yNode(n1) y2 = yNode(n2) y3 = yNode(n3) y4 = yNode(n4) vx1 = v(1, n1) vx2 = v(1, n2) vx3 = v(1, n3) vx4 = v(1, n4) vy1 = v(2, n1) vy2 = v(2, n2) vy3 = v(2, n3) vy4 = v(2, n4) delVX = fhival(s, vx1, vx2) - fhival(s, vx4, vx3) delVY = fhival(s, vy1, vy2) - fhival(s, vy4, vy3) s_dp = s angle = chord(fArg(1, i), s_dp, fArg(2, i)) unitX = COS(angle) unitY = SIN(angle) unitBX = SIN(angle) unitBY = -COS(angle) spread = delVX * unitBX + delVY * unitBY sinist = delVX * unitX + delVY * unitY dExtra = -sinist ! Note conversion of units to mm/year in next line: rLt = dExtra * unitL / unitT opens = spread closes = -spread dip = fhival(s, fDip(1, i), fDip(2, i)) vertic = (ABS(dip - 1.57079632679490D0) < wedge) IF (vertic) THEN relVZ = 0.0D0 probVZ = 1.0D0 ELSE ! Note conversion of units to mm/year in next line: relVZ = closes * ABS(TAN(dip)) * unitL / unitT CALL Likely (delVZ(1, j), delVZ(2, j), relVZ, & ! INTENT(IN) & probVZ) ! INTENT(OUT) END IF CALL Likely (rLat(1, j), rLat(2, j), rLt, & ! INTENT(IN) & probRL) ! INTENT(OUT) shrink = shrink * probRL * probVZ rErr = 0.0 IF ((rLat(1, j) /= 0.0).OR.(rLat(2, j) /= 0.0)) THEN nGSDat = nGSDat + 1 IF ((rLt < rLat(1, j)).AND.(rLat(1, j) /= 0.)) THEN nB = nB + 1 rErr = rLat(1, j) - rLt sumB = sumB + rErr sumB2 = sumB2 + rErr**2 bigB = MAX(bigB, rErr) END IF IF ((rLt > rLat(2, j)).AND.(rLat(2, j) /= 0.)) THEN nB = nB + 1 rErr = rLt - rLat(2, j) sumB = sumB + rErr sumB2 = sumB2 + rErr**2 bigB = MAX(bigB, rErr) END IF END IF zErr = 0.0D0 IF ((delVZ(1, j) /= 0.0D0).OR.(delVZ(2, j) /= 0.0D0)) THEN nGSDat = nGSDat + 1 IF ((relVZ < delVZ(1, j)).AND.(delVZ(1, j) /= 0.0D0)) THEN nB = nB + 1 zErr = delVZ(1, j) - relVZ sumB = sumB + zErr sumB2 = sumB2 + zErr**2 bigB = MAX(bigB, zErr) END IF IF ((relVZ > delVZ(2, j)).AND.(delVZ(2, j) /= 0.0D0)) THEN nB = nB + 1 zErr = relVZ - delVZ(2, j) sumB = sumB + zErr sumB2 = sumB2 + zErr**2 bigB = MAX(bigB, zErr) END IF END IF ! Convert some numbers to '?' or "---" or "agrees" : IF (rLat(1, j) == 0.0D0) THEN c8r1 = ' ?' ELSE WRITE(c8r1, "(F8.3)") rLat(1, j) END IF IF (rLat(2, j) == 0.0D0) THEN c8r2 = ' ?' ELSE WRITE(c8r2, "(F8.3)") rLat(2, j) END IF IF (delVZ(1, j) == 0.0D0) THEN c8v1 = ' ?' ELSE WRITE(c8v1, "(F8.3)") delVZ(1, j) END IF IF (delVZ(2, j) == 0.0D0) THEN c8v2 = ' ?' ELSE WRITE(c8v2, "(F8.3)") delVZ(2, j) END IF IF ((rLat(1, j) == 0.0D0).AND.(rLat(2, j) == 0.0D0)) THEN c8rerr = ' ?' ELSE IF (rErr == 0.0D0) THEN c8rerr = ' agrees' ELSE WRITE (c8rerr, "(F8.3)") rErr END IF IF ((delVZ(1, j) == 0.0D0).AND.(delVZ(2, j) == 0.0D0)) THEN c8verr = ' ?' ELSE IF (zErr == 0.0D0) THEN c8verr = ' agrees' ELSE WRITE (c8verr, "(F8.3)") zErr END IF IF ((rLat(1, j) == 0.0D0).AND.(rLat(2, j) == 0.0D0)) THEN c7rcha = ' ------' ELSE WRITE (c7rcha, "(F7.4)") probRL END IF IF ((delVZ(1, j) == 0.0D0).AND.(delVZ(2, j) == 0.0D0)) THEN c7vcha = ' ------' ELSE WRITE (c7vcha, "(F7.4)") probVZ END IF IF (vertic) THEN WRITE(iUnitT, 1370) fName(j)(1:100), & & c8r1, c8r2, & & rLt, c8rerr, c7rcha 1370 FORMAT(' ', A, 2A8, F8.3, A8, A7, & & ' Throw not used with strike-slip fault.') IF (table_of_offset_rates) THEN WRITE (51, 1370) fName(j)(1:100), c8r1, c8r2, rLt, c8rerr, c7rcha END IF ! table_of_offset_rates (optional log file when running interactively) ELSE WRITE(iUnitT, 1380) fName(j)(1:100), & & c8r1, c8r2, & & rLt, c8rerr, c7rcha, & & c8v1, c8v2, & & relVZ, c8verr, c7vcha 1380 FORMAT(' ', A, 2A8, F8.3, A8, A7, 2A8, F8.3, A8, A7) IF (table_of_offset_rates) THEN WRITE (51, 1380) fName(j)(1:100), c8r1, c8r2, rLt, c8rerr, c7rcha, c8v1, c8v2, relVZ, c8verr, c7vcha END IF ! table_of_offset_rates (optional log file when running interactively) END IF 1390 CONTINUE ! Summary geologic slip-rate statistics for this model: WRITE (iUnitT, 1392) TRIM(title1), TRIM(title2), TRIM(title3) 1392 FORMAT (/' FOR MODEL WITH TITLES:'/ & & ' ', A / ' ', A / ' ', A /) enB = FLOAT(nB) / FLOAT(nGSDat) sumB = sumB / FLOAT(nGSDat) sumB2 = SQRT(sumB2 / FLOAT(nGSDat)) WRITE (iUnitT, 1394) enB, sumB, sumB2, bigB, shrink 1394 FORMAT ( & & /' FRACTION of predictions in error =', F6.3,'.' & & /' AVERAGE error of all predictions = ', F10.3, ' mm/year.' & & /' Root-Mean-Square error of all predictions =', F10.3, ' mm/year.' & & /' LARGEST error of any prediction = ', F10.3, ' mm/year.' & & /' Probability that model input parameters and physical equations are consistent' & & /' with the true geology represented by this dataset is', D10.3, '.') chance = shrink**(1.0D0 / nGSDat) WRITE (iUnitT, 1398) chance 1398 FORMAT( & & ' ',11('====')/ & & ' # Estimate, based on this dataset, of #'/ & & ' # probability that, in a given case, #'/ & & ' # these physical equations and parameter #'/ & & ' # values lead to a predicted slip rate #'/ & & ' # consistent with geologic data (within #'/ & & ' # its uncertainties) is ',F8.5,'. #'/ & & ' ', 11('====')/) IF (table_of_offset_rates) THEN WRITE (51, 1392) TRIM(title1), TRIM(title2), TRIM(title3) WRITE (51, 1394) enB, sumB, sumB2, bigB, shrink WRITE (51, 1398) chance CLOSE (51) END IF ! table_of_offset_rates (optional log file when running interactively) ! Use root-mean-square error (in mm/year) as summary statistic: slpErr = sumB2 ! ****** END GEOLOGIC SLIP-RATE SCORING ******** !------------------------------------------------------------------- ! ****** SPREADING-RATE SCORING SECTION ******** ! Read seafloor-spreading rates on mid-ocean ridges for scoring: 501 CALL GetMOR (iUnitM, iUnitT, mxMOR, & ! INTENT(IN) & MORTag, MORPhi, MORThe, & ! INTENT(OUT) & MORVel, MORSig, numMOR) IF (numMOR <= 0) THEN spread = 0.0D0 GO TO 600 END IF toler = 200.0D0 / 6371.0D0 ! maximum spatial discrepancy, in radians (= km/km) WRITE (iUnitT, 530) 530 FORMAT (// & &' Seafloor-spreading velocities versus model predictions:'// & &' Ridge E.Lon. N.Lat. Velocity (mm/yr) Model.V (mm/yr)', & &' Error (mm/yr) (Sigmas)'/ & &' ----- ------ ------ -------- ----- ------- -----', & &' ----- ----- ------') 531 FORMAT (' ', A5, F7.2, F7.2, ES9.2, ' (', F5.1, ')', ES10.2, & &' (', F5.1, ')', ES9.2, ' (', F5.1,') (', F4.1, ')') sum0s = 0.0D0 sum1 = 0.0D0 sum1s = 0.0D0 sum2 = 0.0D0 sum2s = 0.0D0 sumN = 0.0D0 sumNs = 0.0D0 numDen = numMOR DO 590 i = 1, numMOR ! For each mid-ocean ridge spreading-rate, theLon = oezOPi * MORPhi(i) theLat = 90.0 - oezOPi * MORThe(i) ! Find closest normal fault: rS(1) = SIN(MORThe(i)) * COS(MORPhi(i)) rS(2) = SIN(MORThe(i)) * SIN(MORPhi(i)) rS(3) = COS(MORThe(i)) rMin = 9.99E29 DO 550 j = 1, nFl dip = 0.5D0 * (fDip(1, j) + fDip(2, j)) ! Dip must be greater than 37.5 degrees: IF (ABS(dip - 1.57079632679490D0) < 0.916294D0) THEN ! Dip must be less than 64.0 degrees: IF (ABS(dip - 1.57079632679490D0) > 0.448585D0) THEN ! Positions of end-nodes of fault: rA(1) = SIN(xNode(nodeF(1, j))) * & & COS(yNode(nodeF(1, j))) rA(2) = SIN(xNode(nodeF(1, j))) * & & SIN(yNode(nodeF(1, j))) rA(3) = COS(xNode(nodeF(1, j))) rB(1) = SIN(xNode(nodeF(2, j))) * & & COS(yNode(nodeF(2, j))) rB(2) = SIN(xNode(nodeF(2, j))) * & & SIN(yNode(nodeF(2, j))) rB(3) = COS(xNode(nodeF(2, j))) ! Vectpr along fault from node 1 to 2: vf(1) = rB(1) - rA(1) vf(2) = rB(2) - rA(2) vf(3) = rB(3) - rA(3) fLong = SQRT(vf(1)**2 + vf(2)**2 + vf(3)**2) ! Vector from node 1 to datum: voff(1) = rS(1) - rA(1) voff(2) = rS(2) - rA(2) voff(3) = rS(3) - rA(3) ! Determine internal coordinate (0 ... 1) of datum: frac = (voff(1) * vf(1) + voff(2) * vf(2) + voff(3) * vf(3)) / (fLong**2) IF (frac < 0.0D0) THEN frac = 0.0D0 r = SQRT(voff(1)**2 + voff(2)**2 + voff(3)**2) ELSE IF (frac > 1.0D0) THEN frac = 1.0D0 r = SQRT((rS(1) - rB(1))**2 + & & (rS(2) - rB(2))**2 + & & (rS(3) - rB(3))**2) ELSE CALL DCross (vf, voff, vt) vtl = SQRT(vt(1)**2 + vt(2)**2 + vt(3)**2) r = vtl / fLong END IF IF (r < rMin) THEN rMin = r iEle = j fSave = frac END IF END IF END IF 550 CONTINUE IF (rMin > toler) THEN WRITE (iUnitT, 551) i, MORTag(i), theLon, theLat, toler 551 FORMAT (' #', I5, ' ON RISE ', A5,' AT ',F7.2, & & 'E, ',F6.2,'N IS', & & ' MORE THAN ',F5.3,' RADIANS', & & ' FROM NORMAL FAULT: IGNORED.') numDen = numDen - 1 ELSE ! Compute spreading-rate (iEle, fSave) dvat = v(1, nodeF(4, iEle)) - v(1, nodeF(1, iEle)) dvap = v(2, nodeF(4, iEle)) - v(2, nodeF(1, iEle)) dvbt = v(1, nodeF(3, iEle)) - v(1, nodeF(2, iEle)) dvbp = v(2, nodeF(3, iEle)) - v(2, nodeF(2, iEle)) dvt = dvat + fSave * (dvbt - dvat) dvp = dvap + fSave * (dvbp - dvap) ! Find theta and phi unit-vectors at this point: uP(1) = -rS(2) uP(2) = rS(1) uP(3) = 0.0D0 CALL Unit (uP) ! INTENT(INOUT) CALL DCross (uP, rS, uT) voff(1) = dvt * uT(1) + dvp * uP(1) voff(2) = dvt * uT(2) + dvp * uP(2) voff(3) = dvt * uT(3) + dvp * uP(3) rA(1) = SIN(xNode(nodeF(1, iEle))) * & & COS(yNode(nodeF(1, iEle))) rA(2) = SIN(xNode(nodeF(1, iEle))) * & & SIN(yNode(nodeF(1, iEle))) rA(3) = COS(xNode(nodeF(1, iEle))) rB(1) = SIN(xNode(nodeF(2, iEle))) * & & COS(yNode(nodeF(2, iEle))) rB(2) = SIN(xNode(nodeF(2, iEle))) * & & SIN(yNode(nodeF(2, iEle))) rB(3) = COS(xNode(nodeF(2, iEle))) ! Vector along fault from node 1 to 2: vf(1) = rB(1) - rA(1) vf(2) = rB(2) - rA(2) vf(3) = rB(3) - rA(3) fLong = SQRT(vf(1)**2 + vf(2)**2 + vf(3)**2) CALL DCross (vf, voff, vt) dot = vt(1) * rS(1) + vt(2) * rS(2) + vt(3) * rS(3) spread = dot / fLong bad = ABS(spread - MORVel(i)) badSig = bad / MORSig(i) mvmma = MORVel(i) * secPYr * 1000.0D0 sprmma = spread * secPYr * 1000.0D0 badMma = bad * secPYr * 1000.0D0 WRITE (iUnitT, 531) MORTag(i), theLon, theLat, & & MORVel(i), mvmma, & & spread, sprmma, bad, badMma, badSig IF (badSig > 2.0D0) sum0s = sum0s + 1 sum1 = sum1 + bad sum1s = sum1s + badSig sum2 = sum2 + bad**2 sum2s = sum2s + badSig**2 sumN = MAX(sumN, bad) sumNs = MAX(sumNs, badSig) END IF 590 CONTINUE IF (numDen > 0) THEN sum0s = sum0s / numDen sum1 = sum1 / numDen ELSE sum0s = 0.0D0 sum1 = 0.0D0 END IF sum1mm = sum1 * 1000.0D0 * secPYr IF (numDen > 0) THEN sum1s = sum1s / numDen sum2 = SQRT(sum2 / numDen) ELSE sum1s = 0.0D0 sum2 = 0.0D0 END IF sum2mm = sum2 * secPYr * 1000.0D0 IF (numDen > 0) THEN sum2s = SQRT(sum2s / numDen) ELSE sum2s = 0.0D0 END IF sumNmm = 1000. * secPYr * sumN WRITE (iUnitT, 599) numDen, sum0s, sum1, sum1mm, sum1s, & & sum2, sum2mm, sum2s, sumN, sumNmm, sumNs 599 FORMAT(/' Summary of spreading prediction errors:' & & /' Number of spreading-rates used for scoring: ', I6 & & /' Fraction of predictions wrong by over 2 sigma: ', F5.3 & & /' Mean velocity error: ', ES10.2, ' (', F7.2, ' mm/yr)' & & /' Mean number of sigmas in error: ', F7.2, & & /' RMS velocity error : ', ES10.2, ' (', F7.2, ' mm/yr)' & & /' RMS number of sigmas in error: ', F5.2, & & /' Worst velocity error :', ES10.2, ' (', F7.2, ' mm/yr)' & & /' Worst error, in sigmas: ', F7.2) ! Arbitrarily select an error measure (mean error, in mm/year): spread = sum1mm ! *********** END SPREADING-RATE SCORING ************ !---------------------------------------------------------------- ! *********** BEGIN SEISMICITY SCORING ************ 600 numEqs = 0 time1 = 9999.0D0 * secPYr time2 = 0.0D0 ! Attempt to read catalog file (which may not be present) ! in format produced by my Seismicity.f90 program: WRITE(iUnitT, 601) iUnitE 601 FORMAT(/' ATTEMPTING TO READ SEISMIC CATALOG FILE ON UNIT', I3/) 602 READ(iUnitE, 603, IOSTAT = ios) flag10, & & jYear, month, kDay, & & iHour, minute, iSecon, iTenth, & & eqELon(numEqs + 1), eqNLat(numEqs + 1), & & iDepth, eqMag(numEqs + 1) 603 FORMAT(A10, & & I4, 1X, I2, 1X, I2, 1X, & & I2, 1X, I2, 1X, I2, 1X, I1, 1X, & & F8.3, 1X, F7.3, 1X, & & I3, F6.2) IF (ios /= 0) THEN ! Problem with last READ (e.g., no file, or end-of-file) GO TO 610 ELSE numEqs = numEqs + 1 IF (numEqs > maxEqs) THEN WRITE(iUnitT, 604) 604 FORMAT(/' INCREASE PARAMETER maxEqs TO EQUAL THE NUMBER OF EARTHQUAKES IN INPUT FILE.') CALL Pause() STOP END IF tSec = jYear * secPYr + (month - 1) * 2.6298D6 + (kDay - 1) * 86400.0D0 + & & (iHour - 1) * 3600.0D0 + (minute - 1) * 60.0D0 + iSecon + iTenth / 10.0D0 ! Time1 and time2 are earliest and last earthquake times, ! in seconds. time1 = MIN(time1, tSec) time2 = MAX(time2, tSec) GO TO 602 END IF 610 seismi = 0.0D0 ! (Initializing error measure in case of no catalog) IF (numEqs <= 0) GO TO 700 WRITE(iUnitT, 611) numEqs 611 FORMAT(//' -----------------------------------------------' & & /' Beginning seismicity scoring section;' & & //' Read catalog file, found ', I8, ' earthquakes.') ! Smooth delta-function seismicity with a Gaussian spatial filter, ! and scale to convert to a scalar strain-rate; ! evaluate this strain-rate at each node (for plotting): ! ======================= sigma = 60.0D0 * 1000.0D0 ! ======================= ! Sigma is the smoothing distance, in meters. ! Note that it should allow for epicenter mislocation and ! also for finite source dimensions. ! However, it is not used to achieve a regionally-smooth ! scalar strain-rate field; this will be done below. sigma2 = sigma**2 vP = 8100.0D0 rho = 3300.0D0 ! (SI values appropriate for upper mantle/oceanic lithosphere) vS = vP / SQRT(3.) SIMu = rho * vS**2 ! SIMu is the elastic shear modulus, in Pa: thick = 20000.0D0 ! Thick is the thickness of the seismogenic layer, in m. deltaT = time2 - time1 ! deltaT is the length of the catalog, in seconds. factor = 1.0D0 / (3.14159265358979D0 * SIMu * sigma2 * thick * deltaT) DO 615 i = 1, numNod eDotNC(i) = 0.0D0 uvecN(1, i) = SIN(xNode(i)) * COS(yNode(i)) uvecN(2, i) = SIN(xNode(i)) * SIN(yNode(i)) uvecN(3, i) = COS(xNode(i)) 615 CONTINUE DO 630 j = 1, numEqs IF ((flag10(1:3) == 'ISC').OR.(flag10(1:3) == 'IGN')) THEN ! ISC or IGN catalog, so magnitudes are body-wave ! magnitudes or m_b: eqMB = eqMag(j) eqMS = (eqMB - 2.5D0) / 0.63D0 ! Conversion from m_b to surface-wave magnitude (eqMS) ! from Richter's Elementary Seismology textbook. eqM0 = 10.0D0**(1.5D0 * eqMS + 9.14D0) ! Conversion from M_s to moment (eqM0) from ! Ekstrom & Dziewonski (1998, Nature, vol. 332, pp. 319-323), ! which updates formula of ! Kanamori and Anderson (BSSA, vol. 65, p. 1073), ! in which the constant was +8.5. ! The formala gives eqM0 in Newton-meters (SI). ELSE IF ((flag10(1:10) == 'HarvardCMT').OR.(flag10(1:9) == 'GlobalCMT')) THEN ! Magnitudes are moment-magnitudes already. eqMW = eqMag(j) eqM0 = 10.0D0**(1.5D0 * eqMW + 9.1D0) ! per definition of moment-magnitude in ! Hanks and Kanamori (1979, J. Geophys. Res., pp. 2348-2350). ELSE IF (flag10(1:2) == 'NZ') THEN ! CHANGED BY Z.L.: NEW ZEALAND LOCAL MAGNITUDE ML IS CALCULATED BASED ON ! FOMULA GIVEN BY HAINES IN BSSA,71,275-294,1981. SYSTEMATIC DIFFERENCE ! BETWEEN ML AND BODY MAGNITUDE mb EXISTS AND IS DEPTH-DEPENDENT. SEE AL ! DAVID HARTE AND DAVID VERE-JONES, N.Z.J.GEOL.GEOPHYS.,42,237-253,1999 ! VERY ROUGHLY, mb-ML ~0.2 FOR EVENTS SHALLOWER THAN 40KM BASED ON HISTO ! GIVEN IN HARTE AND VERE-JONES' PAPER. SO WE SIMPLY ADD 0.2 FROM ML TO ! IT INTO mb. THEN CONVERT mb INTO MS AND MW USING ABOVE FORMULA eqMB = eqMag(j) + 0.2D0 eqMS = (eqMB - 2.5D0) / 0.63D0 eqM0 = 10.0D0**(1.5D0 * eqMS + 9.14D0) ELSE WRITE (iUnitT, 617) TRIM(flag10) 617 FORMAT(/' ILLEGAL CATALOG FLAG: ', A) CALL Pause() STOP ' ' END IF equVec(1) = COS(eqNLat(j) * piO180) * & & COS(eqELon(j) * piO180) equVec(2) = COS(eqNLat(j) * piO180) * & & SIN(eqELon(j) * piO180) equVec(3) = SIN(eqNLat(j) * piO180) DO 620 i = 1, numNod dot = equVec(1) * uvecN(1, i) + & & equVec(2) * uvecN(2, i) + & & equVec(3) * uvecN(3, i) CALL DCross (equVec, uvecN(1, i), cUvec) croSiz = SQRT(cUvec(1)**2 + cUvec(2)**2 + cUvec(3)**2) arc = ATAN2(croSiz, dot) r2 = (arc * radius)**2 arg = -r2 / sigma2 IF (arg > -80.0D0) THEN ! (Avoid underflows for distant nodes) eDotNC(i) = eDotNC(i) + factor * eqM0 * EXP(arg) END IF 620 CONTINUE 630 CONTINUE ! Average catalog strain-rate to element centers ! and impose reasonable minimum of 0.1% / 5 Ga = 6.0D-21 /s ! (to prevent correlation test of logarithms from ! being dominated by artificial, uncertain low-end values!) DO 640 i = 1, numEl eDotEC(i) = (eDotNC(nodes(1, i)) + & & eDotNC(nodes(2, i)) + & & eDotNC(nodes(3, i))) / 3.0D0 eDotEC(i) = MAX(eDotEC(i), 6.0D-21) 640 CONTINUE ! END SEISMIC-CATALOG SCALAR STRAIN-RATE FIELD (ROUGH!) ! ------------------------------------------------------------ ! BEGIN FINITE-ELEMENT-MODEL SCALAR STRAIN-RATE FIELD ! Modify velocity solution by averaging together the velocities ! of all nodes that share a common location, thus forcing the ! velocity field to be continuous (and eliminating infinite ! strain-rates in fault elements): DO 642 i = 1, numNod checkN(i) = .FALSE. ! (which means "not yet involved in averaging') 642 CONTINUE DO 670 i = 1, nFl DO 660 j1 = 1, 2 nj1 = nodeF(j1, i) ! (Fault ends are the only places that can have delta.V's) IF (.NOT.checkN(nj1)) THEN list(1) = nj1 checkN(nj1) = .TRUE. ! Begin list of neighbors with paired node: j2 = 5 - j1 nj2 = nodeF(j2, i) list(2) = nj2 checkN(nj2) = .TRUE. nInSum = 2 ! Find shortest fault connected to either one: dx = xNode(nj1) - xNode(nj2) dy = yNode(nj1) - yNode(nj2) IF (dy > 3.14159265358979D0) dy = dy - 6.28318530717959D0 IF (dy < -3.14159265358979D0) dy = dy + 6.28318530717959D0 dy = dy * SIN(xNode(nj1)) short = SQRT(dx**2 + dy**2) DO 646 k = 1, nFl nl1 = nodeF(1, k) nl2 = nodeF(2, k) nl3 = nodeF(3, k) nl4 = nodeF(4, k) IF ((nj1 == nl1).OR.(nj2 == nl1).OR. & & (nj1 == nl2).OR.(nj2 == nl2).OR. & & (nj1 == nl3).OR.(nj2 == nl3).OR. & & (nj1 == nl4).OR.(nj2 == nl4)) THEN dx = xNode(nl1) - xNode(nl2) dy = yNode(nl1) - yNode(nl2) IF (dy > 3.14159265358979D0) dy = dy - 6.28318530717959D0 IF (dy < -3.14159265358979D0) dy = dy + 6.28318530717959D0 dy = dy * SIN(xNode(nl1)) test = SQRT(dx**2 + dy**2) short = MIN(short, test) END IF 646 CONTINUE ! Collect all corner nodes within 10% of this distance: toler = short / 10.0D0 t2 = toler**2 DO 650 k = 1, numNod IF (.NOT.checkN(k)) THEN dx = xNode(nj1) - xNode(k) dy = yNode(nj1) - yNode(k) IF (dy > 3.14159265358979D0) dy = dy - 6.28318530717959D0 IF (dy < -3.14159265358979D0) dy = dy + 6.28318530717959D0 dy = dy * SIN(xNode(nj1)) r2 = dx**2 + dy**2 IF (r2 < t2) THEN nInSum = nInSum + 1 IF (nInSum > mxStar) THEN WRITE(iUnitT, 648) 648 FORMAT(/' INCREASE VALUE' & & ,' OF PARAMETER maxAtP.') CALL Pause() STOP END IF list(nInSum) = k checkN(k) = .TRUE. END IF END IF 650 CONTINUE ! Averaging actually begins here, based on list: xSum = 0.0D0 ySum = 0.0D0 DO 652 k = 1, nInSum xSum = xSum + v(1, list(k)) ySum = ySum + v(2, list(k)) 652 CONTINUE xMean = xSum / nInSum yMean = ySum / nInSum DO 656 k = 1, nInSum v(1, list(k)) = xMean v(2, list(k)) = yMean 656 CONTINUE END IF 660 CONTINUE 670 CONTINUE ! Compute model strain-rates after removing delta.V's: CALL EDot (dxs, dys, & & fpsfer, mxEl, & & mxNode, nodes, numEl, radius, sita, v, & ! INTENT(IN) & eRate) ! INTENT(OUT) ! Convert element-center values to scalars: ! ((e3*-e1*)/2) at integration point 1 ! (Note: * means strain-rate partitioned(?) if E2 /= 0) DO 672 i = 1, numEl exx = eRate(1, 1, i) eyy = eRate(2, 1, i) exy = eRate(3, 1, i) CALL Prince (exx, eyy, exy, & ! INTENT(IN) & e1, e2, u1x, u1y, u2x, u2y) ! INTENT(OUT) ez = -(exx + eyy) IF ((e2 == 0.0D0).AND.(e1 == 0.0D0)) THEN eDotEM(i) = 1.0D-20 ELSE IF ((e2 * ez) > 0.0D0) THEN ! e1 has the unique sign and is partitioned: e1Part = .TRUE. e2Part = .FALSE. eZPart = .FALSE. ELSE IF ((e1 * ez) > 0.0D0) THEN ! e2 has the unique sign and is partitioned: e1Part = .FALSE. e2Part = .TRUE. eZPart = .FALSE. ELSE ! ezz has the unique sign and is partitioned: e1Part = .FALSE. e2Part = .FALSE. eZPart = .TRUE. END IF ! strike-slip rate (e2me1) IF (e1Part) THEN e2me1 = 2.0D0 * ABS(e2) ELSE IF (e2Part) THEN e2me1 = 2.0D0 * ABS(e1) ELSE e2me1 = ABS(e2 - e1) END IF ! thrust-faulting rate (ezme1) IF (e1Part) THEN ezme1 = 2.0D0 * ABS(ez) ELSE IF (eZPart) THEN ezme1 = 2.0D0 * ABS(e1) ELSE ezme1 = ABS(ez - e1) END IF ! normal-faulting rate (e2mez) IF (e2Part) THEN e2mez = 2.0D0 * ABS(ez) ELSE IF (eZPart) THEN e2mez = 2.0D0 * ABS(e2) ELSE e2mez = ABS(e2 - ez) END IF eLarge = MAX(e2me1, ezme1, e2mez) eDotEM(i) = MAX(eLarge / 2.0D0, 1.0D-20) END IF 672 CONTINUE ! Repeatedly smooth scalar strain-rate fields ! by averaging element values onto nodes ! (treating faults as healed), ! and then interpolating from nodes to element centers. ! The same smoothing operator is applied to both ! catalog and finite-element-model scalar strain-rate fields. ! ==================================== nBland = 16 ! ==================================== nBland = MAX(nBland, 1) WRITE(iUnitT, 673) nBland 673 FORMAT(/' Using ', I2, ' smoothing cycles on both strain-rates.') DO 684 nB = 1, nBland ! (A) Extrapolate from elements to nodes (across faults): DO 674 i = 1, numNod eDotNC(i) = 0.0D0 eDotNM(i) = 0.0D0 atnode(i) = 0.0D0 674 CONTINUE DO 676 i = 1, numEl DO 675 j = 1, 3 node = nodes(j, i) eDotNC(node) = eDotNC(node) + eDotEC(i) eDotNM(node) = eDotNM(node) + eDotEM(i) atnode(node) = atnode(node) + 1.0D0 675 CONTINUE 676 CONTINUE DO 678 i = 1, nFl DO 677 j = 1, 2 node1 = nodeF(j, i) node2 = nodeF(5 - j, i) t1 = eDotNC(node1) t2 = eDotNC(node2) eDotNC(node1) = eDotNC(node1) + t2 eDotNC(node2) = eDotNC(node2) + t1 t1 = eDotNM(node1) t2 = eDotNM(node2) eDotNM(node1) = eDotNM(node1) + t2 eDotNM(node2) = eDotNM(node2) + t1 t1 = atnode(node1) t2 = atnode(node2) atnode(node1) = atnode(node1) + t2 atnode(node2) = atnode(node2) + t1 677 CONTINUE 678 CONTINUE DO 679 i = 1, numNod eDotNC(i) = eDotNC(i) / MAX(atnode(i), 1.0D0) eDotNM(i) = eDotNM(i) / MAX(atnode(i), 1.0D0) 679 CONTINUE ! (B) Interpolate from nodes to element centers: DO 680 i = 1, numEl eDotEC(i) = (eDotNC(nodes(1, i)) + & & eDotNC(nodes(2, i)) + & & eDotNC(nodes(3, i))) / 3.0D0 eDotEM(i) = (eDotNM(nodes(1, i)) + & & eDotNM(nodes(2, i)) + & & eDotNM(nodes(3, i))) / 3.0D0 680 CONTINUE 684 CONTINUE ! ----- Output smoothed scalar strain-rate fields ! to two .FEG files: ---------------------- WRITE(c3, "(I3)") nBland WRITE(iUnitT, 686) sigma, iUnitO 686 FORMAT(/' EQ magnitudes (m_b) have been converted to moments,' & & /' the moments distributed with a Gaussian filter of' & & /' characteristic distance ', ES10.3,' m;' & & /' then the scalar strain-rate was diffusively smoothed;' & & /' Use the .feg file written to unit ', I3 & & ,' to plot the results.'/) enctit = 'STRAIN-RATES FROM SEISMIC CATALOG, ' // 'diffusively smoothed (n =' // c3 // ')' CALL PutNet (iUnitO, & & .TRUE., eDotNC, eDotNC, fDip, & & mxEl, mxFEl, mxNode, n1000, & & nFakeN, nFl, nodeF, nodes, & & nRealN, numEl, numNod, offset, & & enctit, eDotNC, xNode, yNode, eDotNC) ! INTENT(IN) WRITE(iUnitT, 688) iUnitZ 688 FORMAT(/' Finite-element model velocities were made continuous' & & /' by healing all active faults;' & & /' then the scalar strain-rate was diffusively smoothed;' & & /' Use the .feg file written to unit ', I3 & & ,' to plot the results.'/) enctit = 'STRAIN-RATES FROM HEALED-FAULT FE MODEL, ' // 'diffusively smoothed (n =' // c3 // ')' CALL PutNet (iUnitZ, & & .TRUE., eDotNM, eDotNM, fDip, & & mxEl, mxFEl, mxNode, n1000, & & nFakeN, nFl, nodeF, nodes, & & nRealN, numEl, numNod, offset, & & enctit, eDotNM, xNode, yNode, eDotNM) ! INTENT(IN) ! ----------------------------------------------- ! Compute correlation coefficient between: ! *log10(smoothed scalar strain-rate from seismic catalog) ! *log10(smoothed scalar strain-rate from fault-healed model) ! Convert to common logarithms: DO 690 i = 1, numEl eDotEC(i) = LOG10(MAX(eDotEC(i), 1.0D-30)) eDotEM(i) = LOG10(MAX(eDotEM(i), 1.0D-30)) 690 CONTINUE ! Find averages of log10's of scalar strain-rates; catalog & model: aveC = 0.0D0 aveM = 0.0D0 sum = 0.0D0 DO 693 i = 1, numEl aveC = aveC + area(i) * eDotEC(i) aveM = aveM + area(i) * eDotEM(i) sum = sum + area(i) 693 CONTINUE aveC = aveC / sum aveM = aveM / sum varC = 0.0D0 varM = 0.0D0 crossP = 0.0D0 DO 694 i = 1, numEl varC = varC + area(i) * (eDotEC(i) - aveC)**2 varM = varM + area(i) * (eDotEM(i) - aveM)**2 crossP = crossP + area(i) * (eDotEC(i) - aveC) * (eDotEM(i) - aveM) 694 CONTINUE correl = crossP / SQRT(varC * varM) WRITE(iUnitT, 695) correl 695 FORMAT(/' Correlation coefficient between:' & & /' * log10(smoothed scalar strain-rate from seismic catalog)' & & /' * log10(smoothed scalar strain-rate from fault-headed model)' & & /' is: ', F7.4) ! Arbitrarily select an error measure (correlation coefficient) seismi = correl ! *********** END SEISMICITY SCORING ************ !---------------------------------------------------------------- ! ****** UPPER-MANTLE SEISMIC ANISOTROPY SCORING SECTION ******** ! N.B. The basic assumption is that upper-mantle seismic anisotropy ! (expressed as SKS fast-polarization {"phi"} azimuths and ! splitting times in s} is due to simple shear in the ! asthenosphere, and thus is comparable to directions of ! basal shear tractions inferred from a torque-report file ! (but only for slab-less plates, with slab_q(iPlate) = .FALSE.). ! Bear in mind that anisotropy could also be due to several other causes: ! (1) Most-extensional direction of Neogene finite strain in lithosphere. ! (2) Most-extensional direction of ancient finite strain in lithosphere. ! (3) Azimuth of aligned magma-filled dikes (or water-filled cracks). ! and in these cases the data would NOT be comparable as assumed! ! READ fast-polarization azimuths ("phi") and splitting times ! of SKS waves (e.g., compilation of Matt Fouch). 700 CALL GetSKS (iUnitK, iUnitT, mxSKS, & & numSKS, SKS_tag, SKS_theta, SKS_phi, & ! INTENT(IN) & SKS_argument, SKS_delay) ! INTENT(OUT) IF (numSKS <= 0) THEN anisotropy = 0.0D0 ! summary error measure GO TO 900 END IF WRITE (iUnitT, 701) numSKS, iUnitK 701 FORMAT (/' ',I6,' Fast SKS azimuths and delay times were read from unit ', I2) !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - IF (DIs_Initialization_Needed()) THEN subdivision = 4 WRITE (*, "(/' Initializing uniform global grid at subdivision ', & & I1, ' for area-weighting of data...')") subdivision CALL DInitialize_Weighting (subdivision) END IF WRITE (*, "(' Computing area-weights for SKS azimuth/delay data...')") IF (ALLOCATED(weights)) DEALLOCATE ( weights ) ALLOCATE ( weights(numSKS) ) CALL DPerform_Weighting (number_of_data = numSKS, & & theta_radians = SKS_theta, & & phi_radians = SKS_phi, & ! INTENT(IN) & weights = weights) ! INTENT(OUT) !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ! Read plate outlines, for assigning sites to the correct plate: WRITE (*, 2) iUnitF 2 FORMAT (/' ATTEMPTING TO READ OUTLINES OF PLATES FROM UNIT', I3/) CALL GetPBx (iUnitF, iUnitT, names, nPBnd, nPlate, & ! INTENT(IN) & nBoundaryPoints, pLat, pLon) ! INTENT(OUT) !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - WRITE (*, "(' Assigning SKS azimuth/delay data to plates...')") DO 720 i = 1, numSKS CALL ThetaPhi_2_pLate (iUnitT, & & nPBnd, nBoundaryPoints, & & nPlate, pLat, pLon, & & SKS_theta(i), SKS_phi(i), & ! INTENT(IN) & iPlate) ! INTENT(OUT) SKS_iPlate(i) = iPlate 720 CONTINUE !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - WRITE (*, "(' Finding basal shear azimuths to compare with SKS azimuths...')") CALL Tractor(iunitQ, iUnitT, nPlate, numSKS, & & slab_q, SKS_phi, SKS_theta, SKS_iPlate, & ! INTENT(IN) & basal_shear_tractions) ! INTENT(OUT) !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ! Prepare to compute azimuths of fastest-horizontal-extension principal strain-rate ! axes, as alternate predictor of fast-polarization azimuths (phi) of SKS data: ! Compute the (theta, phi) and then (x, y, z) of integration points: DO 750 i = 1, numEl DO 725 k = 1, 3 CartVs(1, k) = COS(yNode(nodes(k, i))) * & & SIN(xNode(nodes(k, i))) CartVs(2, k) = SIN(yNode(nodes(k, i))) * & & SIN(xNode(nodes(k, i))) CartVs(3, k) = COS(xNode(nodes(k, i))) 725 CONTINUE DO 740 m = 1, 7 DO 730 j = 1, 3 tempV(j) = 0.0 DO 728 k = 1, 3 tempV(j) = tempV(j) + CartVs(j, k) * points(k, m) 728 CONTINUE 730 CONTINUE CALL Unit (tempV) ! INTENT(INOUT) IF (ABS(tempV(3)) <= 0.5D0) THEN xIP(m, i) = ACOS(tempV(3)) ELSE equPar = SQRT(tempV(1)**2 + tempV(2)**2) xIP(m, i) = ATan2F(equPar, tempV(3)) END IF yIP(m, i) = ATan2F(tempV(2), tempV(1)) CartR(1, m, i) = SIN(xIP(m, i)) * COS(yIP(m, i)) CartR(2, m, i) = SIN(xIP(m, i)) * SIN(yIP(m, i)) CartR(3, m, i) = COS(xIP(m, i)) 740 CONTINUE 750 CONTINUE ! Strain rates at integration points: CALL EDot (dxs, dys, & & fpsfer, mxEl, & & mxNode, nodes, numEl, radius, sita, v, & ! INTENT(IN) & eRate) ! INTENT(OUT) ! Note: Integration point must be within one-half-element- ! width of datum, or datum is not used. ! Here, the typical element side is determined by ! averaging together all element sides, ! and expressed in radians: sumSid = 0.0D0 DO 756 i = 1, numEl DO 753 k = 1, 3 kp1 = 1 + MOD(k, 3) n1 = nodes(k, i) n2 = nodes(kp1, i) x1 = xNode(n1) x2 = xNode(n2) y1 = yNode(n1) y2 = yNode(n2) rA(1) = SIN(x1) * COS(y1) rA(2) = SIN(x1) * SIN(y1) rA(3) = COS(x1) rB(1) = SIN(x2) * COS(y2) rB(2) = SIN(x2) * SIN(y2) rB(3) = COS(x2) dot = rA(1) * rB(1) + rA(2) * rB(2) + rA(3) * rB(3) CALL DCross (rA, rB, cUvec) croSiz = SQRT(cUvec(1)**2 + cUvec(2)**2 + cUvec(3)**2) side = ATan2F(croSiz, dot) sumSid = sumSid + side 753 CONTINUE 756 CONTINUE toler = sumSid / (6.0D0 * numEl) !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - WRITE (iUnitT, 760) 760 FORMAT (//' Fast SKS-polarization azimuths versus model predictions (in slab-less plates):'/ & &/' Datum E.Long. N.Latt. Plate Delay,s Azimuth Model.Az1,2 Error Area-weight Slabless plate?' & &/' ----- ------- ------- ----- ------- ------- ----------- ----- ----------- ---------------') 761 FORMAT ( ' ', A5, 1X, F7.2, 1X, F7.2, 4X, A2, 1X, F7.2, 1X, F7.0, 1X, F5.0, ',', F5.0, 1X, F5.0, 1X, F11.4, 1X, ' .TRUE. ') 762 FORMAT ( ' ', A5, 1X, F7.2, 1X, F7.2, 4X, A2, 1X, F7.2, 1X, F7.0, 1X, ' ??, ??', 1X, ' --', 9X, '---', 1X, ' .FALSE. ') !: Initialize accumulators: nSKS_bad = 0 ! SKS data misfit by more than 20 degrees (L0 numerator) nSKS_used = 0 ! Count of data in slabless plates (L0 denominator) sum1 = 0.0D0 ! Numerator for L1 (mean size of error) sum1d = 0.0D0 ! Denominator for L1 sum2 = 0.0D0 ! Numerator for L2 (RMS size of error) sum2d = 0.0D0 ! Denominator for L2 DO 790 i = 1, numSKS ! Next SKS datum: theLon = oezOPi * SKS_phi(i) IF (theLon > 180.0D0) theLon = theLon - 360.0D0 IF (theLon < -180.0D0) theLon = theLon + 360.0D0 theLat = 90.0D0 - oezOPi * SKS_theta(i) !Cartesian unit vector: rS(1) = SIN(SKS_theta(i)) * COS(SKS_phi(i)) rS(2) = SIN(SKS_theta(i)) * SIN(SKS_phi(i)) rS(3) = COS(SKS_theta(i)) !Datum: Azimuth of fast polarization of SKS, in degrees: azimuth = 180.0D0 - oezOPi * SKS_argument(i) IF (azimuth < 0.0D0) azimuth = azimuth + 180.0D0 IF (azimuth < 0.0D0) azimuth = azimuth + 180.0D0 IF (azimuth >= 180.0D0) azimuth = azimuth - 180.0D0 IF (azimuth >= 180.0D0) azimuth = azimuth - 180.0D0 IF (slab_q(SKS_iPlate(i))) THEN !No scoring; just echo the datum, without providing any model predictions... WRITE (iUnitT, 762) SKS_tag(i), theLon, theLat, names(SKS_iPlate(i)), SKS_delay(i), azimuth ELSE ! slabless-plate; proceed with scoring! !Predictor1: Azimuth of simple shear in asthenosphere, in degrees: bAzim = 180.0D0 - oezOPi * ATan2F(basal_shear_tractions(2, i), basal_shear_tractions(1, i)) ! (y = phi = E, x = theta = S) IF (bAzim < 0.0D0) bAzim = bAzim + 180.0D0 IF (bAzim < 0.0D0) bAzim = bAzim + 180.0D0 IF (bAzim >= 180.0D0) bAzim = bAzim - 180.0D0 IF (bAzim >= 180.0D0) bAzim = bAzim - 180.0D0 !Predictor2: Azimuth of e3 (fastest-extension) horizontal principal strain-rate in lithosphere, in degrees: !(a) Find closest integration point in the grid... r2min = 9.99D29 DO 770 m = 1, 7 DO 765 j = 1, numEl r2 = (rS(1) - CartR(1, m, j))**2 + & & (rS(2) - CartR(2, m, j))**2 + & & (rS(3) - CartR(3, m, j))**2 IF (r2 < r2min) THEN iEle = j mIP = m r2min = r2 END IF 765 CONTINUE 770 CONTINUE rMin = SQRT(r2min) IF (rMin > toler) THEN WRITE (iUnitT, 780) strtag(i), theLon, theLat, toler 780 FORMAT (' ', A5, 1X, F7.2, 1X, F7.2, & & ' SKS DATUM IS MORE THAN ', & & F6.4, ' RADIANS FROM NEAREST' & & ,' INTEGRATION POINT: IGNORED.') e3Azim = bAzim ! so prediction errors will be the same, and bAzim error won't be underbid. ELSE !(b) Get e_1H (fastest-horizontal-shortening) direction: e11h = eRate(1, mIP, iEle) e22h = eRate(2, mIP, iEle) e12h = eRate(3, mIP, iEle) CALL Prince (e11h, e22h, e12h, & ! INTENT(IN) & e1h, e2h, u1x, u1y, u2x, u2y) ! INTENT(OUT) err = -(e11h + e22h) !(or, = -(e1h + e2h) ) e1 = MIN(e1h, e2h, err) e3 = MAX(e1h, e2h, err) IF ((err > e1).AND.(err < e3)) THEN e2 = err ELSE IF ((e1h > e1).AND.(e1h < e3)) THEN e2 = e1h ELSE e2 = e2h END IF pAzim = 180.0D0 - oezOPi * ATan2F(u1y, u1x) IF (pAzim < 0.0D0) pAzim = pAzim + 180.0D0 IF (pAzim < 0.0D0) pAzim = pAzim + 180.0D0 IF (pAzim >= 180.0D0) pAzim = pAzim - 180.0D0 IF (pAzim >= 180.0D0) pAzim = pAzim - 180.0D0 e3Azim = pAzim + 90.0D0 END IF IF (e3Azim < 0.0D0) e3Azim = e3Azim + 180.0D0 IF (e3Azim < 0.0D0) e3Azim = e3Azim + 180.0D0 IF (e3Azim >= 180.0D0) e3Azim = e3Azim - 180.0D0 IF (e3Azim >= 180.0D0) e3Azim = e3Azim - 180.0D0 !Prediction error of Predictor1 (simple shear in asthenosphere): simple_bad = ABS(bAzim - azimuth) IF (simple_bad > 90.0D0) simple_bad = ABS(simple_bad - 180.0D0) IF (simple_bad > 90.0D0) simple_bad = ABS(simple_bad - 180.0D0) IF (simple_bad > 90.0D0) simple_bad = ABS(simple_bad - 180.0D0) IF (simple_bad > 90.0D0) simple_bad = ABS(simple_bad - 180.0D0) !Prediction error of Predictor2 (pure-shear in lithosphere): pure_bad = ABS(e3Azim - azimuth) IF (pure_bad > 90.0D0) pure_bad = ABS(pure_bad - 180.0D0) IF (pure_bad > 90.0D0) pure_bad = ABS(pure_bad - 180.0D0) IF (pure_bad > 90.0D0) pure_bad = ABS(pure_bad - 180.0D0) IF (pure_bad > 90.0D0) pure_bad = ABS(pure_bad - 180.0D0) !Select smaller of two prediction errors: bad = MIN(simple_bad, pure_bad) IF (bad > 20.0D0) nSKS_bad = nSKS_bad + 1 ! SKS data misfit by more than 20 degrees (L0 numerator) nSKS_used = nSKS_used + 1 ! Count of data in slabless plates (L0 denominator) sum1 = sum1 + bad * SKS_delay(i) * weights(i) ! Numerator for L1 (mean size of error) sum1d = sum1d + SKS_delay(i) * weights(i) ! Denominator for L1 sum2 = sum2 + bad**2 * SKS_delay(i) * weights(i) ! Numerator for L2 (RMS size of error) sum2d = sum2d + SKS_delay(i) * weights(i) ! Denominator for L2 WRITE (iUnitT, 761) SKS_tag(i), theLon, theLat, names(SKS_iPlate(i)), SKS_delay(i), azimuth, bAzim, e3Azim, bad, weights(i) END IF ! slab_q(SKS_iPlate(i)), or NOT <== proceed with scoring 790 CONTINUE IF (nSKS_used > 0) THEN perBad = (100.0D0 * nSKS_bad) / (1.0D0 * nSKS_used) ELSE perBad = 0.0D0 END IF IF (sum1d > 0.0D0) THEN sum1 = sum1 / sum1d ELSE sum1 = 0.0D0 END IF IF (sum2d > 0.0D0) THEN sum2 = SQRT(sum2 / sum2d) ELSE sum2 = 0.0D0 END IF WRITE (iUnitT, 796) nSKS_used, perBad, sum1, sum2 796 FORMAT (/' SUMMARY OF SEISMIC ANISOTROPY ERR0RS:' & & /' Number of SKS data in slabless plates, used for scoring: ', I6, '.' & & /' Percentage of directions mis-predicted by more than 20 degrees: ', F6.2 & & /' Area- & delay-weighted MEAN error: ', F5.2, ' degrees.' & & /' Area- & delay-weighted RMS error: ', F5.2, ' degrees.') ! Arbitrarily choose measure of error (mean error, in degrees): anisotropy = sum1 DEALLOCATE (weights) ! ********* END ANISOTROPY SCORING SECTION ********** !---------------------------------------------------------------- ! ****** SUMMARY OF SCORES ******** 900 WRITE (iUnitT, 901) IF (geodes /= 0.0D0) WRITE (iUnitT, 910) geodes IF (spread /= 0.0D0) WRITE (iUnitT, 920) spread IF (stress /= 0.0D0) WRITE (iUnitT, 930) stress IF (slpErr /= 0.0D0) WRITE (iUnitT, 940) slpErr IF (seismi /= 0.0D0) WRITE (iUnitT, 950) seismi IF (anisotropy /= 0.0D0) WRITE (iUnitT, 970) anisotropy 901 FORMAT (' ------------------------------------------------' & & //' REVIEW OF SUMMARY MEASURES OF ERROR OR QUALITY:') 910 FORMAT (' GEODETIC VELOCITY (ERROR): ', F7.2) 920 FORMAT (' SEAFLOOR SPREADING RATE (ERROR): ', F7.2) 930 FORMAT (' STRESS DIRECTION (ERROR): ', F7.2) 940 FORMAT (' FAULT SLIP RATE (ERROR): ', F7.2) 950 FORMAT (' (SMOOTHED) SEISMICITY CORRELATION (QUALITY): ', F7.4) 970 FORMAT (' SEISMIC ANISOTROPY (ERROR): ', F7.2) WRITE (iUnitT, 999) TRIM(title1), TRIM(title2), TRIM(title3) 999 FORMAT (/' FOR MODEL:' & & /' ', A & & /' ', A & & /' ', A/) CALL Pause() ! <=== NOTE: Check in source code for SUBROUTINE Pause() to see whether ! user-response ( [Enter] key) is required, or skipped? ! If program does not wait for user [Enter], it will ! terminate before you can read its summary statistics! ! (This is desirable behavior when program is run in batch mode, ! which is usually the case in the middle of serious modeling experiments.) STOP CONTAINS SUBROUTINE Adjust (data, eLon, iUnitT, mxAdj, n, nLat, & & radius, weights, & ! INTENT(IN) & predic, & ! INTENT(INOUT) & eLonP, nLatP, rate) ! INTENT(OUT) ! Compares predicted (horizontal) velocities "predic" ! (whose first component is the Southward or +Theta velocity, ! and whose second component is the Eastward or +Phi velocity) ! with geodetic horizontal velocities "data" (in the same format) ! which have assigned (dimensionless) "weights" of mean value 1.0, ! at a set of points i=1, ... ,"n" <= "mxAdj" ! (defined by "eLon(i)" and "nLat(i)", in degrees East & North) ! and finds the rotation correction (to add) which minimizes the ! RMS of the magnitude of the velocity differences. ! It then adds this correction to "predic", and reports ! the pole in user-friendly coordinates "eLonP" and "nLatP" ! (in degrees East and North) ! and gives the rotation rate "rate" in radians/second. USE DSphere ! in Peter Bird's file DSphere.f90 (just for radians_per_degree, degrees_per_radian) USE MKL95_PRECISION ! Intel's Math Kernel Library USE MKL95_LAPACK ! Intel's Math Kernel Library, LAPACK (Linear Analysis Package) portion IMPLICIT NONE INTEGER, INTENT(IN) :: iUnitT, mxAdj, n REAL*8, DIMENSION(2, mxAdj), INTENT(IN) :: data ! components (Theta, Phi) = (S, E) REAL*8, DIMENSION(mxAdj), INTENT(IN) :: eLon, nLat, weights REAL*8, INTENT(IN) :: radius REAL*8, DIMENSION(2, mxAdj), INTENT(INOUT) :: predic ! components (Theta, Phi) = (S, E) REAL*8, INTENT(OUT) :: eLonP, nLatP, rate CHARACTER*1 :: uplo INTEGER :: i, info, lwork INTEGER, DIMENSION(:), ALLOCATABLE :: ipiv REAL*8 :: DeltaVPhi, DeltaVTheta, dPhi, dTheta, dx, dy, dz, & & equat, & & pPhi, pTheta, px, py, pz, & & rx, ry, rz, & & w REAL*8, DIMENSION(3) :: DeltaVXyz, dXyz, pXyz, rMeters, uPhi, uTheta, ur REAL*8, DIMENSION(3, 3) :: coef REAL*8, DIMENSION(3, 1) :: right ! N.B. Array "right" requires 2 subscripts, although only the first column is used. ! This is an abritrary rule imposed by MKL routine dsysv, which only accepts 2-subscript arrays as inputs. REAL*8, DIMENSION(3) :: solut REAL*8, DIMENSION(:), ALLOCATABLE :: work IF (n > mxAdj) THEN WRITE (iUnitT, 1) mxAdj, n 1 FORMAT (//' ERR0R: PARAMETER mxAdj (NOW ', I8, ')' & & /' MUST BE INCREASED TO AT LEAST ', I8) CALL Pause() STOP END IF coef = 0.0D0 ! all (1:3, 1:3) components right = 0.0D0 ! all (1:3, 1:1) components DO i = 1, n ! build up coefficient matrix and right-hand-side of linear system by summing over benchmarks: CALL DLonLat_2_Uvec(eLon(i), nLat(i), ur) CALL DLocal_Theta(ur, uTheta) ! N.B. This will crash if any geodetic benchmark is exactly on the N or S pole. CALL DLocal_Phi (ur, uPhi) ! N.B. This will crash if any geodetic benchmark is exactly on the N or S pole. rx = radius * ur(1) ry = radius * ur(2) rz = radius * ur(3) dTheta = data(1, i) dPhi = data(2, i) dXyz(1:3) = dTheta * uTheta(1:3) + dPhi * uPhi(1:3) dx = dXyz(1) dy = dXyz(2) dz = dXyz(3) pTheta = predic(1, i) pPhi = predic(2, i) pXyz(1:3) = pTheta * uTheta(1:3) + pPhi * uPhi(1:3) px = pXyz(1) py = pXyz(2) pz = pXyz(3) w = weights(i) !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - coef(1, 1) = coef(1, 1) + w * ( 2.0D0 * rz**2 + 2.0D0 * ry**2) ! coefficient of Omega_x in equation #1 (minimization w.r.t. Omega_x) coef(1, 2) = coef(1, 2) + w * (-2.0D0 * ry * rx) ! coefficient of Omega_y in equation #1 (minimization w.r.t. Omega_x) coef(1, 3) = coef(1, 3) + w * (-2.0D0 * rz * rx) ! coefficient of Omega_z in equation #1 (minimization w.r.t. Omega_x) right(1, 1) = right(1, 1) - w * (-2.0D0 * rz * (py - dy) + 2.0D0 * ry * (pz - dz)) ! right-hand side of equation #1 (minimization w.r.t. Omega_x) !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - coef(2, 1) = coef(2, 1) + w * (-2.0D0 * rx * ry) ! coefficient of Omega_x in equation #2 (minimization w.r.t. Omega_y) coef(2, 2) = coef(2, 2) + w * ( 2.0D0 * rz**2 + 2.0D0 * rx**2) ! coefficient of Omega_y in equation #2 (minimization w.r.t. Omega_y) coef(2, 3) = coef(2, 3) + w * (-2.0D0 * rz * ry) ! coefficient of Omega_z in equation #2 (minimization w.r.t. Omega_y) right(2, 1) = right(2, 1) - w * ( 2.0D0 * rz * (px - dx) - 2.0D0 * rx * (pz - dz)) ! right-hand side of equation #2 (minimization w.r.t. Omega_y) !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - coef(3, 1) = coef(3, 1) + w * (-2.0D0 * rx * rz) ! coefficient of Omega_x in equation #3 (minimization w.r.t. Omega_z) coef(3, 2) = coef(3, 2) + w * (-2.0D0 * ry * rz) ! coefficient of Omega_y in equation #3 (minimization w.r.t. Omega_z) coef(3, 3) = coef(3, 3) + w * ( 2.0D0 * ry**2 + 2.0D0 * rx**2) ! coefficient of Omega_z in equation #3 (minimization w.r.t. Omega_z) right(3, 1) = right(3, 1) - w * (-2.0D0 * ry * (px - dx) + 2.0D0 * rx * (py - dy)) ! right-hand side of equation #3 (minimization w.r.t. Omega_z) !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - END DO !Option to write-out this linear system before attempting a solution: WRITE (iUnitT, *) WRITE (iUnitT, "(' Linear system to be solved in subprogram Adjust:')") DO i = 1, 3 IF (i == 2) THEN WRITE (iUnitT, "(/ ' ', 3ES11.3, ' X Omega_', I1, ' = ', ES11.3)") (coef(i, j), j = 1, 3), i, right(i, 1) ELSE WRITE (iUnitT, "(/ ' ', 3ES11.3, ' Omega_', I1, ' ', ES11.3)") (coef(i, j), j = 1, 3), i, right(i, 1) END IF END DO WRITE (iUnitT, *) !Using a routine from Intel's Math Kernel Library, Linear Analysis Package (LAPACK) portion: !Using dsysv for this SYMMETRIC INDEFINITE linear system: uplo = 'U' ! problem is stated in the Upper Triangle (plus diagonal) of the coefficient matrix. ALLOCATE ( ipiv(3) ) !lwork = -1 ! First time: "lwork = -1" signals a workspace-size-query to dsysv; answer to be placed in work(1). Answer was "3.0000000D0". lwork = 3 ! Based on test described above. ALLOCATE ( work(lwork) ) !call dsysv(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info ) ! <=== This is the sample Fortran77 CALL in the manual. CALL dsysv(uplo, 3, 1, coef, 3, ipiv, right, 3, work, lwork, info) ! N.B. I'm using the Fortran77 CALL because the F95 CALL is buggy. IF (info /= 0) THEN WRITE (iUnitT, "(' Error when Adjust calls dsysv of MKL_LAPACK: info = ', I12)") info CALL Pause() STOP END IF DEALLOCATE ( work ) DEALLOCATE ( ipiv ) solut(1:3) = right(1:3, 1) ! unpack the solution. Note that only the first column has been used. !Characterize the Euler pole for the rotation that will be added: rate = SQRT(solut(1)**2 + solut(2)**2 + solut(3)**2) IF (rate > 0.0D0) THEN equat = SQRT(solut(1)**2 + solut(2)**2) eLonP = degrees_per_radian * ATAN2(solut(2), solut(1)) nLatP = degrees_per_radian * ATAN2(solut(3), equat) !Add this rotation to all velocities in predic: DO i = 1, n CALL DLonLat_2_Uvec(eLon(i), nLat(i), ur) CALL DLocal_Theta(ur, uTheta) ! N.B. This will crash if any geodetic benchmark is exactly on the N or S pole. CALL DLocal_Phi (ur, uPhi) ! N.B. This will crash if any geodetic benchmark is exactly on the N or S pole. rMeters(1:3) = radius * ur(1:3) CALL DCross(solut, rMeters, DeltaVXyz) DeltaVTheta = DDot(DeltaVXyz, uTheta) DeltaVPhi = DDot(DeltaVXyz, uPhi) predic(1, i) = predic(1, i) + DeltaVTheta predic(2, i) = predic(2, i) + DeltaVPhi END DO ELSE eLonP = 0.0D0 nLatP = 0.0D0 END IF END SUBROUTINE Adjust REAL*8 FUNCTION ATanPV (y, x) ! Returns principal value of inverse tangent of (y, x). ! N.B. Order of these two arguments is NOT an error; ! this order corresponds to ("motivator", "restrainer"). IMPLICIT NONE REAL*8, INTENT(IN) :: y, x REAL*8 :: xx, yy IF (x >= 0.0D0) THEN xx = x yy = y ELSE xx = -x yy = -y END IF IF ((yy == 0.0D0).OR.(xx == 0.0D0)) THEN IF (yy == 0.0D0) ATanPV = 0.0D0 IF (xx == 0.0D0) ATanPV = 1.57079632679490D0 ELSE ATanPV = ATAN2(yy, xx) END IF END FUNCTION ATanPV REAL*8 FUNCTION Chord (angle1, s, angle2) ! Returns an angle obtained by interpolation between angle1 ! and angle2. The interpolation method is NOT sensitive to ! possible cycle shifts (of 2*n*Pi) between angle1 and angle2. ! Unit vectors are constructed for angle1 and angle2, and a ! linear chord is drawn between their tips. ! S is the internal coordinate along the chord; ! it is dimensionless, with value 0.0D0 at angle1 and 1.0D0 at ! angle2. (The user may input s values outside this range ! to get results outside the (smaller) angle between angle1 and ! angle2, if desired.) The angle returned is that from the ! origin to this chord point. ! This algorithm should work equally well for angles measured ! either clockwise or counterslockwise from any reference, as ! long as the usage is consistent. ! Both the input angles and the result "chord" are in radians. IMPLICIT NONE REAL*8, INTENT(IN) :: angle1, s, angle2 REAL*8 :: c1, c2 REAL*8, DIMENSION(2) :: uvec1, uvec2, vecs uvec1(1) = COS(angle1) uvec1(2) = SIN(angle1) uvec2(1) = COS(angle2) uvec2(2) = SIN(angle2) vecs(1) = (1.0D0 - s) * uvec1(1) + s * uvec2(1) vecs(2) = (1.0D0 - s) * uvec1(2) + s * uvec2(2) c1 = vecs(1) c2 = vecs(2) Chord = ATan2F(c2, c1) END FUNCTION Chord SUBROUTINE Deriv (iUnitT, mxEl, mxNode, & ! input & nodes, numEl, & & radius, xNode, yNode, & & area, detJ, & ! output & dXS, dYS, dXSP, dYSP, fPSfer, sita) ! Sets up 6 vector nodal functions (fPSfer) of each spherical ! triangle finite element, at each of its 7 integration points. ! Calculates dXS and dYS, the Theta-derivitive and Phi-derivitive ! of each of these 6 vector nodal functions. ! Also computes "area", the areas of the plane triangles. ! Also computes "detJ", the local ratio of areas on the sphere ! to areas on the plane triangles. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - INTEGER, INTENT(IN) :: iUnitT, mxEl, mxNode, nodes, numEl ! input REAL*8, INTENT(IN) :: radius, xNode, yNode ! input REAL*8, INTENT(OUT) :: area, detJ, dXS, dYS, dXSP, dYSP, fPSfer, sita ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - DOUBLE PRECISION points DOUBLE PRECISION fff, skkc, skke, sncsne, snccse, csccse, cscsne DOUBLE PRECISION xa, xb, xc, ya, yb, yc, za, zb, zc, xyzp INTEGER i, j, m REAL*8 a, areaP, b, c, cka, cosm, cscs, cse, cssn, & & dd, dd1, dd2, dd3, ddpn, dpdc, dpde, & & pfq, phaij, phi, pnx, pny, pnz, pp, rn, rr1, rr2, rr3, & & sitaj, sitami, snc, sne, theta, tx, ty, & & x21, x31, y21, y31, z21, z31 DIMENSION xNode(mxNode), yNode(mxNode), nodes(3, mxEl), area(mxEl) DIMENSION detJ(7, mxEl) DIMENSION dXS(2, 2, 3, 7, mxEl), dYS(2, 2, 3, 7, mxEl) DIMENSION dXSP(3, 7, mxEl), dYSP(3, 7, mxEl), points(3, 7) DIMENSION phi(3), theta(3), skkc(3), skke(3), fff(3), & & sita(7, mxEl), fPSfer(2, 2, 3, 7, mxEl) COMMON / S1S2S3 / points DO 900 i = 1, numEl DO 100 j = 1, 3 theta(j) = xNode(nodes(j, i)) phi(j) = yNode(nodes(j, i)) 100 CONTINUE x21 = SIN(theta(2)) * COS(phi(2)) - SIN(theta(1)) * COS(phi(1)) x31 = SIN(theta(3)) * COS(phi(3)) - SIN(theta(1)) * COS(phi(1)) y21 = SIN(theta(2)) * SIN(phi(2)) - SIN(theta(1)) * SIN(phi(1)) y31 = SIN(theta(3)) * SIN(phi(3)) - SIN(theta(1)) * SIN(phi(1)) z21 = COS(theta(2)) - COS(theta(1)) z31 = COS(theta(3)) - COS(theta(1)) a = y21 * z31 - y31 * z21 b = z21 * x31 - z31 * x21 c = x21 * y31 - x31 * y21 areaP = SQRT(a * a + b * b + c * c) area(i) = radius * radius * (0.5D0 * areaP) pnx = a / areaP pny = b / areaP pnz = c / areaP dd1 = SIN(theta(1)) * COS(phi(1)) * pnx dd2 = SIN(theta(1)) * SIN(phi(1)) * pny dd3 = COS(theta(1)) * pnz dd = dd1 + dd2 + dd3 ! This part is to test if Kong's method and Bird's method give the same ! results for the derivitive: xa = SIN(theta(1)) * COS(phi(1)) xb = SIN(theta(2)) * COS(phi(2)) xc = SIN(theta(3)) * COS(phi(3)) ya = SIN(theta(1)) * SIN(phi(1)) yb = SIN(theta(2)) * SIN(phi(2)) yc = SIN(theta(3)) * SIN(phi(3)) za = COS(theta(1)) zb = COS(theta(2)) zc = COS(theta(3)) cka = (yb * zc - zb * yc) * xa + (zb * xc - xb * zc) * ya + (xb * yc - yb * xc) * za DO 800 m = 1, 7 snccse = 0.0D0 sncsne = 0.0D0 cosm = 0.0D0 DO 200 j = 1, 3 snccse = snccse + points(j, m) * SIN(theta(j)) * COS(phi(j)) sncsne = sncsne + points(j, m) * SIN(theta(j)) * SIN(phi(j)) cosm = cosm + points(j, m) * COS(theta(j)) 200 CONTINUE xyzp = SQRT(snccse * snccse + sncsne * sncsne + cosm * cosm) snccse = snccse / xyzp sncsne = sncsne / xyzp cosm = cosm / xyzp sitaj = ACOS(cosm) ty = sncsne tx = snccse phaij = ATan2F(ty, tx) csccse = COS(sitaj) * COS(phaij) cscsne = COS(sitaj) * SIN(phaij) ! Bird's method: fff(1) = ((yb * zc - zb * yc) * snccse + (zb * xc - xb * zc) * sncsne & & + (xb * yc - yb * xc) * cosm) / cka fff(2) = ((yc * za - zc * ya) * snccse + (zc * xa - xc * za) * sncsne & & + (xc * ya - yc * xa) * cosm) / cka fff(3) = ((ya * zb - za * yb) * snccse + (za * xb - xa * zb) * sncsne & & + (xa * yb - ya * xb) * cosm) / cka skkc(1) = ((yb * zc - zb * yc) * csccse & & + (zb * xc - xb * zc) * cscsne & & - (xb * yc - yb * xc) * SIN(sitaj)) / cka skkc(2) = ((yc * za - zc * ya) * csccse & & + (zc * xa - xc * za) * cscsne & & - (xc * ya - yc * xa) * SIN(sitaj)) / cka skkc(3) = ((ya * zb - za * yb) * csccse & & + (za * xb - xa * zb) * cscsne & & - (xa * yb - ya * xb) * SIN(sitaj)) / cka skke(1) = (-(yb * zc - zb * yc) * sncsne & & + (zb * xc - xb * zc) * snccse) / cka skke(2) = (-(yc * za - zc * ya) * sncsne & & + (zc * xa - xc * za) * snccse) / cka skke(3) = (-(ya * zb - za * yb) * sncsne & & + (za * xb - xa * zb) * snccse) / cka sita(m, i) = sitaj rr1 = SIN(sitaj) * COS(phaij) rr2 = SIN(sitaj) * SIN(phaij) rr3 = COS(sitaj) rn = rr1 * pnx + rr2 * pny + rr3 * pnz pp = dd / rn dpdc = (COS(sitaj) * COS(phaij) * pnx + COS(sitaj) * SIN(phaij) * pny & & - SIN(sitaj) * pnz) dpde = (-SIN(sitaj) * SIN(phaij) * pnx + & & SIN(sitaj) * COS(phaij) * pny) ddpn = pp / rn dpdc = -ddpn * dpdc dpde = -ddpn * dpde IF(sita(m, i) <= 0.0D0.OR.sita(m, i) >= 3.14159265358979D0) THEN sitami = sita(m, i) * 57.2957795130823D0 WRITE(iUnitT, 220) m, i, sitami 220 FORMAT(' COLATITUDE OF INTEGRATION POINT',I5, & & ' OF ELEMENT', & & I5,' IS OUT RANGE', & & D14.4) CALL Pause() STOP END IF DO 500 j = 1, 3 dXSP(j, m, i) = dpdc * fff(j) + pp * skkc(j) dYSP(j, m, i) = dpde * fff(j) + pp * skke(j) cscs = COS(theta(j)) * COS(phi(j)) cssn = COS(theta(j)) * SIN(phi(j)) snc = SIN(theta(j)) sne = SIN(phi(j)) cse = COS(phi(j)) fPSfer(1, 1, j, m, i) = cscs * csccse + cssn * cscsne & & + snc * SIN(sitaj) fPSfer(2, 1, j, m, i) = -sne * csccse + cse * cscsne fPSfer(1, 2, j, m, i) = -cscs * SIN(phaij) + cssn * COS(phaij) fPSfer(2, 2, j, m, i) = sne * SIN(phaij) + cse * COS(phaij) dXS(1, 1, j, m, i) = (-cscs * snccse - cssn * sncsne & & + snc * COS(sitaj)) * fff(j) & & + fPSfer(1, 1, j, m, i) * skkc(j) dXS(2, 1, j, m, i) = (sne * snccse - cse * sncsne) * fff(j) & & + fPSfer(2, 1, j, m, i) * skkc(j) dYS(1, 1, j, m, i) = (-cscs * cscsne + cssn * csccse) * fff(j) & & + fPSfer(1, 1, j, m, i) * skke(j) dYS(2, 1, j, m, i) = (sne * cscsne + cse * csccse) * fff(j) & & + fPSfer(2, 1, j, m, i) * skke(j) dXS(1, 2, j, m, i) = fPSfer(1, 2, j, m, i) * skkc(j) dXS(2, 2, j, m, i) = fPSfer(2, 2, j, m, i) * skkc(j) dYS(1, 2, j, m, i) = (-cscs * COS(phaij) - cssn * SIN(phaij)) & & * fff(j) & & + fPSfer(1, 2, j, m, i) * skke(j) dYS(2, 2, j, m, i) = (sne * COS(phaij) - cse * SIN(phaij)) & & * fff(j) & & + fPSfer(2, 2, j, m, i) * skke(j) fPSfer(1, 1, j, m, i) = fPSfer(1, 1, j, m, i) * fff(j) fPSfer(2, 1, j, m, i) = fPSfer(2, 1, j, m, i) * fff(j) fPSfer(1, 2, j, m, i) = fPSfer(1, 2, j, m, i) * fff(j) fPSfer(2, 2, j, m, i) = fPSfer(2, 2, j, m, i) * fff(j) 500 CONTINUE pfq = fff(1) + fff(2) + fff(3) ! orphan statement, left over from some test? (pfq does not seem to be used.) detJ(m, i) = rn**3 / (dd * dd) 800 CONTINUE 900 CONTINUE END SUBROUTINE Deriv SUBROUTINE Deriv3 (iUnitT, mxEl, mxNode, & & nodes, numEl, & & radius, xNode, yNode, & ! INTENT(IN) & area, detJ, & ! INTENT(OUT) & dxs, dys, dxsp, dysp, fpsfer, sita) ! Sets up 6 vector nodal functions (fpsfer) of each spherical ! triangle finite element, at each of its 7 integration points. ! Calculates dxs and dys, the theta-derivitive and phi-derivitive ! of each of these 6 vector nodel functions. ! Also computes area, the areas of the plane triangles. ! Also computes detJ, the local ratio of areas on the sphere ! to areas on the plane triangles. USE DSphere ! in Peter Bird's file DSphere.f90 IMPLICIT NONE INTEGER, INTENT(IN) :: iUnitT, mxEl, mxNode, numEl INTEGER, DIMENSION(3, mxEl), INTENT(IN) :: nodes REAL*8, INTENT(IN) :: radius REAL*8, DIMENSION(mxNode), INTENT(IN) :: xNode, yNode REAL*8, DIMENSION(mxEl), INTENT(OUT) :: area REAL*8, DIMENSION(7, mxEl), INTENT(OUT) :: detJ REAL*8, DIMENSION(2, 2, 3, 7, mxEl), INTENT(OUT) :: dxs, dys REAL*8, DIMENSION(3, 7, mxEl), INTENT(OUT) :: dxsp, dysp REAL*8, DIMENSION(2, 2, 3, 7, mxEl), INTENT(OUT) :: fpsfer REAL*8, DIMENSION(7, mxEl), INTENT(OUT) :: sita INTEGER :: i, j, m REAL*8 :: a, areap, b, & & c, cka, cosm, csccse, cscs, cscsne, cse, cssn, & & dd, dd1, dd2, dd3, ddpn, dpdc, dpde, & & fff(3), & & phaij, pfq, phi(3), pnx, pny, pnz, pp, & & rn, rr1, rr2, rr3, & & sitaj, sitami, skkc(3), skke(3), snc, sne, sncsne, snccse, & & theta(3), tx, ty, & & x21, x31, xa, xb, xc, xyzp, & & y21, y31, ya, yb, yc, & & z21, z31, za, zb, zc REAL*8, DIMENSION(3, 7) :: points COMMON / S1S2S3 / points DO 900 i = 1, numEl DO 100 j = 1, 3 theta(j) = xNode(nodes(j, i)) phi(j) = yNode(nodes(j, i)) 100 CONTINUE x21 = SIN(theta(2)) * COS(phi(2)) - SIN(theta(1)) * COS(phi(1)) x31 = SIN(theta(3)) * COS(phi(3)) - SIN(theta(1)) * COS(phi(1)) y21 = SIN(theta(2)) * SIN(phi(2)) - SIN(theta(1)) * SIN(phi(1)) y31 = SIN(theta(3)) * SIN(phi(3)) - SIN(theta(1)) * SIN(phi(1)) z21 = COS(theta(2)) - COS(theta(1)) z31 = COS(theta(3)) - COS(theta(1)) a = y21 * z31 - y31 * z21 b = z21 * x31 - z31 * x21 c = x21 * y31 - x31 * y21 areap = SQRT(a * a + b * b + c * c) area(i) = radius * radius * (0.5 * areap) pnx = a / areap pny = b / areap pnz = c / areap dd1 = SIN(theta(1)) * COS(phi(1)) * pnx dd2 = SIN(theta(1)) * SIN(phi(1)) * pny dd3 = COS(theta(1)) * pnz dd = dd1 + dd2 + dd3 ! This part is to test if Kong's method and Bird's method are same for ! calculating derivative: xa = SIN(theta(1)) * COS(phi(1)) xb = SIN(theta(2)) * COS(phi(2)) xc = SIN(theta(3)) * COS(phi(3)) ya = SIN(theta(1)) * SIN(phi(1)) yb = SIN(theta(2)) * SIN(phi(2)) yc = SIN(theta(3)) * SIN(phi(3)) za = COS(theta(1)) zb = COS(theta(2)) zc = COS(theta(3)) cka = (yb * zc - zb * yc) * xa + (zb * xc - xb * zc) * ya + (xb * yc - yb * xc) * za DO 800 m = 1, 7 snccse = 0.0D0 sncsne = 0.0D0 cosm = 0.0D0 DO 200 j = 1, 3 snccse = snccse + points(j, m) * SIN(theta(j)) * COS(phi(j)) sncsne = sncsne + points(j, m) * SIN(theta(j)) * SIN(phi(j)) cosm = cosm + points(j, m) * COS(theta(j)) 200 CONTINUE xyzp = SQRT(snccse * snccse + sncsne * sncsne + cosm * cosm) snccse = snccse / xyzp sncsne = sncsne / xyzp cosm = cosm / xyzp sitaj = ACOS(cosm) ty = sncsne tx = snccse phaij = ATan2F(ty, tx) csccse = COS(sitaj) * COS(phaij) cscsne = COS(sitaj) * SIN(phaij) ! Bird's method: fff(1) = ((yb * zc - zb * yc) * snccse + (zb * xc - xb * zc) * sncsne & & + (xb * yc - yb * xc) * cosm) / cka fff(2) = ((yc * za - zc * ya) * snccse + (zc * xa - xc * za) * sncsne & & + (xc * ya - yc * xa) * cosm) / cka fff(3) = ((ya * zb - za * yb) * snccse + (za * xb - xa * zb) * sncsne & & + (xa * yb - ya * xb) * cosm) / cka skkc(1) = ((yb * zc - zb * yc) * csccse & & + (zb * xc - xb * zc) * cscsne & & - (xb * yc - yb * xc) * SIN(sitaj)) / cka skkc(2) = ((yc * za - zc * ya) * csccse & & + (zc * xa - xc * za) * cscsne & & - (xc * ya - yc * xa) * SIN(sitaj)) / cka skkc(3) = ((ya * zb - za * yb) * csccse & & + (za * xb - xa * zb) * cscsne & & - (xa * yb - ya * xb) * SIN(sitaj)) / cka skke(1) = (-(yb * zc - zb * yc) * sncsne & & + (zb * xc - xb * zc) * snccse) / cka skke(2) = (-(yc * za - zc * ya) * sncsne & & + (zc * xa - xc * za) * snccse) / cka skke(3) = (-(ya * zb - za * yb) * sncsne & & + (za * xb - xa * zb) * snccse) / cka sita(m, i) = sitaj rr1 = SIN(sitaj) * COS(phaij) rr2 = SIN(sitaj) * SIN(phaij) rr3 = COS(sitaj) rn = rr1 * pnx + rr2 * pny + rr3 * pnz pp = dd / rn dpdc = (COS(sitaj) * COS(phaij) * pnx + COS(sitaj) * SIN(phaij) * pny & & - SIN(sitaj) * pnz) dpde = (-SIN(sitaj) * SIN(phaij) * pnx + & & SIN(sitaj) * COS(phaij) * pny) ddpn = pp / rn dpdc = -ddpn * dpdc dpde = -ddpn * dpde IF ((sita(m, i) <= 0.0D0) .OR. (sita(m, i) >= 3.14159265358979D0)) THEN sitami = sita(m, i) * degrees_per_radian WRITE(iUnitT, 220) m, i, sitami 220 FORMAT(' COLATITUDE OF INTEGRATION POINT', I5, ' OF ELEMENT', I5, ' IS OUT RANGE', E14.4) CALL Pause() STOP END IF DO 500 j = 1, 3 dxsp(j, m, i) = dpdc * fff(j) + pp * skkc(j) dysp(j, m, i) = dpde * fff(j) + pp * skke(j) cscs = COS(theta(j)) * COS(phi(j)) cssn = COS(theta(j)) * SIN(phi(j)) snc = SIN(theta(j)) sne = SIN(phi(j)) cse = COS(phi(j)) fpsfer(1, 1, j, m, i) = cscs * csccse + cssn * cscsne + snc * SIN(sitaj) fpsfer(2, 1, j, m, i) = -sne * csccse + cse * cscsne fpsfer(1, 2, j, m, i) = -cscs * SIN(phaij) + cssn * COS(phaij) fpsfer(2, 2, j, m, i) = sne * SIN(phaij) + cse * COS(phaij) dxs(1, 1, j, m, i) = (-cscs * snccse - cssn * sncsne & & + snc * COS(sitaj)) * fff(j) & & + fpsfer(1, 1, j, m, i) * skkc(j) dxs(2, 1, j, m, i) = (sne * snccse - cse * sncsne) * fff(j) & & + fpsfer(2, 1, j, m, i) * skkc(j) dys(1, 1, j, m, i) = (-cscs * cscsne + cssn * csccse) * fff(j) & & + fpsfer(1, 1, j, m, i) * skke(j) dys(2, 1, j, m, i) = (sne * cscsne + cse * csccse) * fff(j) & & + fpsfer(2, 1, j, m, i) * skke(j) dxs(1, 2, j, m, i) = fpsfer(1, 2, j, m, i) * skkc(j) dxs(2, 2, j, m, i) = fpsfer(2, 2, j, m, i) * skkc(j) dys(1, 2, j, m, i) = (-cscs * COS(phaij) - cssn * SIN(phaij)) & & * fff(j) & & + fpsfer(1, 2, j, m, i) * skke(j) dys(2, 2, j, m, i) = (sne * COS(phaij) - cse * SIN(phaij)) & & * fff(j) & & + fpsfer(2, 2, j, m, i) * skke(j) fpsfer(1, 1, j, m, i) = fpsfer(1, 1, j, m, i) * fff(j) fpsfer(2, 1, j, m, i) = fpsfer(2, 1, j, m, i) * fff(j) fpsfer(1, 2, j, m, i) = fpsfer(1, 2, j, m, i) * fff(j) fpsfer(2, 2, j, m, i) = fpsfer(2, 2, j, m, i) * fff(j) 500 CONTINUE pfq = fff(1) + fff(2) + fff(3) detJ(m, i) = rn**3 / (dd * dd) 800 CONTINUE 900 CONTINUE END SUBROUTINE Deriv3 SUBROUTINE EDot (dXS, dYS, & & fPSfer, mxEl, & & mxNode, nodes, numEl, radius, sita, v, & ! INTENT(IN) & eRate) ! INTENT(OUT) ! Compute strain rate components EDot_xx, EDot_yy, and ! EDot_xy (tensor form; equal to ! (1/2) * ((dVx/dY)+(dVy/dX)) ! at the integration points of triangular continuum elements. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - REAL*8, INTENT(IN) :: dXS, dYS, fPSfer ! input INTEGER, INTENT(IN) :: mxEl, mxNode, nodes, numEl ! input REAL*8, INTENT(IN) :: radius, sita ! input DOUBLE PRECISION, INTENT(IN) :: v ! input REAL*8, INTENT(OUT) :: eRate ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - DOUBLE PRECISION points COMMON / S1S2S3 / points ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - INTEGER i, j, m, node REAL*8 dy11, dy21, dy12, dy22, exx, exy, eyy, fp11, fp21, fp12, fp22, vx, vy DIMENSION dXS(2, 2, 3, 7, mxEl), dYS(2, 2, 3, 7, mxEl), & & eRate(3, 7, mxEl), & & fPSfer(2, 2, 3, 7, mxEl), & & nodes(3, mxEl), points(3, 7), & & sita(7, mxEl), v(2, mxNode) ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - DO 1000 m = 1, 7 DO 900 i = 1, numEl exx = 0.0D0 eyy = 0.0D0 exy = 0.0D0 DO 800 j = 1, 3 node = nodes(j, i) vx = v(1, node) vy = v(2, node) dy11 = dYS(1, 1, j, m, i) / SIN(sita(m, i)) dy21 = dYS(2, 1, j, m, i) / SIN(sita(m, i)) dy12 = dYS(1, 2, j, m, i) / SIN(sita(m, i)) dy22 = dYS(2, 2, j, m, i) / SIN(sita(m, i)) fp11 = fPSfer(1, 1, j, m, i) / TAN(sita(m, i)) fp21 = fPSfer(2, 1, j, m, i) / TAN(sita(m, i)) fp12 = fPSfer(1, 2, j, m, i) / TAN(sita(m, i)) fp22 = fPSfer(2, 2, j, m, i) / TAN(sita(m, i)) exx = exx + vx * dXS(1, 1, j, m, i) + vy * dXS(2, 1, j, m, i) eyy = eyy + vx * dy12 + vy * dy22 + vx * fp11 + vy * fp21 exy = exy + vx * dy11 + vy * dy21 & & + vx * dXS(1, 2, j, m, i) + vy * dXS(2, 2, j, m, i) & & - vx * fp12 - vy * fp22 800 CONTINUE eRate(1, m, i) = exx / radius eRate(2, m, i) = eyy / radius eRate(3, m, i) = 0.50D0 * exy / radius 900 CONTINUE 1000 CONTINUE RETURN END SUBROUTINE EDot SUBROUTINE Extract_LRi (longer_line, & ! input & LRi, shorter_line) ! output ! New routine added for Shells_v5.0+ to support multiple !"Lithospheric Rheology" (abbreviated as "LR") integer codes, ! in any line of the input .feg file which defines an element !(either a triangular continuum element, or a ! linear fault element). ! CHARACTER*80, INTENT(IN) :: longer_line is the whole ! element-definition line from the .feg file. ! INTEGER, INTENT(OUT) :: LRi is the rheologic code ! (or 0, if no such code was found). ! CHARACTER*80, INTENT(OUT) :: shorter_line has the " LRi" portion removed (if any), ! so it can be interpreted by the same code as in Shells_v4.1-. IMPLICIT NONE CHARACTER*80, INTENT(IN) :: longer_line INTEGER, INTENT(OUT) :: LRi CHARACTER*80, INTENT(OUT) :: shorter_line CHARACTER*80 :: string INTEGER :: longer_length, LR_start_byte longer_length = LEN_TRIM(longer_line) LR_start_byte = INDEX(longer_line, "LR") IF (LR_start_byte > 0) THEN ! the "LR" flag was found IF (longer_length > (LR_start_byte + 1)) THEN ! some byte(s) follow the "LR" string = longer_line((LR_start_byte + 2):longer_length) READ (string, *) LRi shorter_line = longer_line(1:(LR_start_byte - 1)) ELSE ! "LR" is present, but nothing follows it; infer 0. LRi = 0 shorter_line = longer_line(1:(LR_start_byte - 1)) END IF ELSE ! no "LR" flag is present LRi = 0 shorter_line = longer_line END IF END SUBROUTINE Extract_LRi SUBROUTINE FAngls (phi, theta, & ! INTENT(IN) & fAngle) ! INTENT(OUT) ! Calculate the arguments (angles counterclockwise from +Theta) ! at both ends of an arc of a great circle. Results in radians. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - REAL*8, INTENT(IN) :: phi, theta ! input REAL*8, INTENT(OUT) :: fAngle ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - DOUBLE PRECISION fPoint COMMON / SFault / fPoint ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - REAL*8 a1, a2, a3, b1, b2, b3, dg180, dx, dxx, dy, dyy, dz, & & p1, p2, phai, s, s1, s2, s3, sita, xval, xx, yy, zz DIMENSION fAngle(2), fPoint(7), phi(2), theta(2) dg180 = 3.14159265358979D0 a1 = SIN(theta(1)) * COS(phi(1)) a2 = SIN(theta(1)) * SIN(phi(1)) a3 = COS(theta(1)) b1 = SIN(theta(2)) * COS(phi(2)) b2 = SIN(theta(2)) * SIN(phi(2)) b3 = COS(theta(2)) s = 0.99D0 xx = s * a1 + (1.0D0 - s) * b1 yy = s * a2 + (1.0D0 - s) * b2 zz = s * a3 + (1.0D0 - s) * b3 xval = SQRT(xx * xx + yy * yy + zz * zz) xx = xx / xval yy = yy / xval zz = zz / xval dx = xx - a1 dy = yy - a2 dz = zz - a3 sita = theta(1) phai = phi(1) s1 = COS(sita) * COS(phai) s2 = COS(sita) * SIN(phai) s3 = -SIN(sita) p1 = -SIN(phai) p2 = COS(phai) dxx = dx * s1 + dy * s2 + dz * s3 dyy = dx * p1 + dy * p2 fAngle(1) = ATan2F(dyy, dxx) s = 0.01D0 xx = s * a1 + (1.0D0 - s) * b1 yy = s * a2 + (1.0D0 - s) * b2 zz = s * a3 + (1.0D0 - s) * b3 xval = SQRT(xx * xx + yy * yy + zz * zz) xx = xx / xval yy = yy / xval zz = zz / xval dx = b1 - xx dy = b2 - yy dz = b3 - zz sita = ACOS(zz) phai = ATan2F(yy, xx) IF(phai < 0.0) phai = 2.0D0 * dg180 + phai s1 = COS(sita) * COS(phai) s2 = COS(sita) * SIN(phai) s3 = -SIN(sita) p1 = -SIN(phai) p2 = COS(phai) dxx = dx * s1 + dy * s2 + dz * s3 dyy = dx * p1 + dy * p2 fAngle(2) = ATan2F(dyy, dxx) RETURN END SUBROUTINE FAngls REAL*8 FUNCTION FltLen (phi1, phi2, radius, theta1, theta2) ! INTENT(IN) ! Calculates length of great circle segment between ! point (theta1, phi1) and point (theta2, phi2), ! in physical length units (radians * radius). ! Note that theta is colatitude (from North pole), ! and phi is East longitude; both in radians. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - REAL*8, INTENT(IN) :: phi1, phi2, radius, theta1, theta2 ! inputs ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - DOUBLE PRECISION ab ab = SIN(theta1) * SIN(theta2) * COS(phi1) * COS(phi2) + & & SIN(theta1) * SIN(theta2) * SIN(phi1) * SIN(phi2) + & & COS(theta1) * COS(theta2) ab = ACOS(ab) FltLen = ab * radius RETURN END FUNCTION FltLen SUBROUTINE FindIt (mxEl, mxNode, nodes, numEl, xNode, yNode, & & theta, phi, & ! INTENT(IN) & iEle, s1, s2, s3) ! INTENT(OUT) ! Determine which element (iEle) contains the surface point ! (theta, phi), and also its internal triangular coordinates ! (s1, s2, s3) in the plate triangle; s1 + s2 + s3 = 1.0D0 ! If the point is not in any element, iEle = 0 is returned. IMPLICIT NONE INTEGER, INTENT(IN) :: mxEl, mxNode, numEl INTEGER, DIMENSION(3, mxEl), INTENT(IN) :: nodes REAL*8, DIMENSION(mxNode), INTENT(IN) :: xNode, yNode REAL*8 :: theta, phi INTEGER, INTENT(OUT) :: iEle REAL*8, INTENT(OUT) :: s1, s2, s3 INTEGER :: i, k REAL*8 :: dot, dotp, eps, growth, & & perp(3), pole(3), rA(3), rB(3), rn(3, 3), rP(3), rS(3), & & s1num, s1den, s2num, s2den DATA eps / 1.0D-4 / ! Cartesian (x, y, z) radius vector to point on sphere: rS(1) = SIN(theta) * COS(phi) rS(2) = SIN(theta) * SIN(phi) rS(3) = COS(theta) DO 100 i = 1, numEl DO 10 k = 1, 3 ! Find Cartesian radius vector for each node: rn(1, k) = SIN(xNode(nodes(k, i))) * COS(yNode(nodes(k, i))) rn(2, k) = SIN(xNode(nodes(k, i))) * SIN(yNode(nodes(k, i))) rn(3, k) = COS(xNode(nodes(k, i))) 10 CONTINUE ! Find normal vector of plane element: rA(1) = rn(1, 2) - rn(1, 1) rA(2) = rn(2, 2) - rn(2, 1) rA(3) = rn(3, 2) - rn(3, 1) rB(1) = rn(1, 3) - rn(1, 2) rB(2) = rn(2, 3) - rn(2, 2) rB(3) = rn(3, 3) - rn(3, 2) CALL DCross (rA, rB, perp) CALL Unit (perp) ! Shorten/lengthen radius vector so it points to plane element: dot = rS(1) * perp(1) + rS(2) * perp(2) + rS(3) * perp(3) IF (dot <= 0.0D0) GO TO 100 dotp = rn(1, 1) * perp(1) + rn(2, 1) * perp(2) + rn(3, 1) * perp(3) growth = dotp / dot rP(1) = rS(1) * growth rP(2) = rS(2) * growth rP(3) = rS(3) * growth ! Compute and test s1, then s2, then s3... CALL DCross (rn(1:3, 2), rn(1:3, 3), pole) s1num = rP(1) * pole(1) + rP(2) * pole(2) + rP(3) * pole(3) s1den = rn(1, 1) * pole(1) + rn(2, 1) * pole(2) + rn(3, 1) * pole(3) IF (s1den /= 0.0D0) THEN s1 = s1num / s1den ELSE GO TO 100 END IF IF ((s1 >= -eps).AND.(s1 <= (1.0D0 + eps))) THEN CALL DCross (rn(1:3, 3), rn(1:3, 1), pole) s2num = rP(1) * pole(1) + rP(2) * pole(2) + rP(3) * pole(3) s2den = rn(1, 2) * pole(1) + rn(2, 2) * pole(2) + rn(3, 2) * pole(3) IF (s2den /= 0.0D0) THEN s2 = s2num / s2den ELSE GO TO 100 END IF IF ((s2 >= -eps).AND.(s2 <= (1.0D0 + eps))) THEN s3 = 1.00D0 - s1 - s2 IF ((s3 >= -eps).AND.(s3 <= (1.0D0 + eps))) THEN iEle = i RETURN END IF END IF END IF 100 CONTINUE iEle = 0 END SUBROUTINE FindIt SUBROUTINE FNodal (mxFEl, mxNode, nFl, nodeF, & ! input & xNode, yNode, & & fPFlt) ! output ! Calculates vector nodal functions at all integration points ! on all arc-of-great-circle fault elements. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - INTEGER, INTENT(IN) :: mxFEl, mxNode, nFl, nodeF ! input REAL*8, INTENT(IN) :: xNode, yNode ! input REAL*8, INTENT(OUT) :: fPFlt ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - DOUBLE PRECISION fPhi COMMON / FPhis / fPhi DIMENSION fPhi(4, 7) ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - INTEGER i, j, k, l, m, n1, n2 REAL*8 fpp, phi, theta DIMENSION fPFlt(2, 2, 2, 7, mxFEl), fpp(2, 2, 2, 7), & & nodeF(4, mxFEl), phi(2), theta(2), & & xNode(mxNode), yNode(mxNode) DO 900 i = 1, nFl n1 = nodeF(1, i) n2 = nodeF(2, i) theta(1) = xNode(n1) theta(2) = xNode(n2) phi(1) = yNode(n1) phi(2) = yNode(n2) CALL SNodal (phi, theta, & ! input & fpp) ! output DO 800 m = 1, 7 DO 500 j = 1, 2 DO 400 k = 1, 2 DO 300 l = 1, 2 fPFlt(l, k, j, m, i) = fpp(l, k, j, m) 300 CONTINUE ! l = 1:2 400 CONTINUE ! k = 1:s 500 CONTINUE ! j = 1:2 800 CONTINUE ! m = 1:7 900 CONTINUE ! i = 1:nFl RETURN END SUBROUTINE FNodal SUBROUTINE GetGEO (iUnitB, iUnitT, iUnitY, iUnitC, mxGeo, & & pltGeo, pltCos, pltInt, & ! INTENT(IN) & geoTag, geoPhi, geoThe, & ! INTENT(OUT) & geoVel, geoSig, geoAzi, & & gpsFMT, numGeo) ! Reads relative horizontal velocity vectors of benchmarks, ! determined by geodesy. ! The input file format expected is: ! Peter Bird's .gps format of 2002.08.01: ! 3 lines of text/format/headers, ! then one line per benchmark, giving ! coordinates, E & N velocity components in mm/a, ! E & N standard deviations of velocity components in mm/a, ! correlation between E and N components of velocity, ! reference frame, and benchmark identification. ! This subprogram will return: ! geoPhi = East longitude (radians) ! geoThe = co-latitude (radians) ! geoVel = velocity magnitude, m/s ! geoAzi = azimuth of velocity, in radians counterclockwise from South. ! geoSig = standard deviation of velocity (assumes circular uncertainty ellipse, and no correlation of errors) ! geoTag = benchmark name USE DSphere ! in Peter Bird's file DSphere.f90 IMPLICIT NONE INTEGER, INTENT(IN) :: iUnitB, iUnitT, iUnitY, iUnitC, mxGeo LOGICAL, INTENT(IN) :: pltGeo, pltCos, pltInt ! Create(?) new .gps-formatted output files for plotting: ! *(Shells model error) and/or ! *(coseismic part of Shells model velocities) and/or ! *(interseismic part of Shells model velocities) ! with -FiniteMap-. CHARACTER*20, DIMENSION(mxGeo), INTENT(OUT) :: geoTag REAL*8, DIMENSION(mxGeo), INTENT(OUT) :: geoPhi, geoThe REAL*8, DIMENSION(mxGeo), INTENT(OUT) :: geoVel, geoSig, geoAzi CHARACTER*132, INTENT(OUT) :: gpsFMT INTEGER, INTENT(OUT) :: numGeo CHARACTER*15 :: frame CHARACTER*20 :: tag CHARACTER*132 :: header INTEGER :: ios REAL*8 :: eLonDg, nLatDg, vEMmpa, vNMmpa, vESigm, vNSigm, correl REAL*8 :: secPYr secPYr = 365.25D0 * 24.0D0 * 60.0D0 * 60.0D0 ! If first READ fails (no file connected), return numGeo = 0 WRITE(iUnitT, 1) iUnitB 1 FORMAT(/' ATTEMPTING TO READ GEODETIC BENCHMARK VELOCITIES (.gps file) FROM UNIT', I3/) numGeo = 0 !Create (optional) output .gps files NOW, so that header lines can be copied into them: IF (pltGeo) OPEN(UNIT = iUnitY, FILE = "ERRORS.gps") ! Create new file, with unconditional OPEN. IF (pltCos) OPEN(UNIT = iUnitC, FILE = "COSEISMIC.gps") ! Create new file, with unconditional OPEN. IF (pltInt) OPEN(UNIT = iUnitI, FILE = "INTERSEISMIC.gps") ! Create new file, with unconditional OPEN. ! Read, but do not use, the first line of text (file mnemonic). READ (iUnitB, * , IOSTAT = ios) IF (ios /= 0) RETURN IF (pltGeo) WRITE(iUnitY, "('ERROR in geodetic velocity predictions')") IF (pltCos) WRITE(iUnitC, "('COSEISMIC part of model geodetic velocities')") IF (pltInt) WRITE(iUnitI, "('INTERSEISMIC part of model geodetic velocities')") ! Read and store the FORMAT in line 2 READ (iUnitB, "(A)", IOSTAT = ios) gpsFMT IF (pltGeo) WRITE(iUnitY, "(A)") TRIM(gpsFMT) IF (pltCos) WRITE(iUnitC, "(A)") TRIM(gpsFMT) IF (pltInt) WRITE(iUnitI, "(A)") TRIM(gpsFMT) ! Read the column headers in line 3 READ (iUnitB, "(A)", IOSTAT = ios) header IF (ios /= 0) RETURN IF (pltGeo) WRITE(iUnitY, "(A)") TRIM(header) IF (pltCos) WRITE(iUnitC, "(A)") TRIM(header) IF (pltInt) WRITE(iUnitI, "(A)") TRIM(header) 10 READ (iUnitB, gpsFMT, END = 100) eLonDg, nLatDg, vEMmpa, vNMmpa, & & vESigm, vNSigm, correl, frame, tag numGeo = numGeo + 1 geoPhi(numGeo) = eLonDg * radians_per_degree geoThe(numGeo) = (90.0D0 - nLatDg) * radians_per_degree geoVel(numGeo) = SQRT(vEMmpa**2 + vNMmpa**2) / (1000.0D0 * secPYr) geoAzi(numGeo) = ATan2F(vEMmpa, -vNMmpa) geoSig(numGeo) = SQRT(vESigm**2 + vNSigm**2) / (1000.0D0 * secPYr) IF (geoSig(numGeo) <= 0.0D0) THEN WRITE (iUnitT, 19) numGeo 19 FORMAT (/' ERROR: Standard deviation of velocity' & & /' is <= 0.0 for benchmark ', I5 & & /' All values must be positive.') CALL Pause () STOP END IF geoTag(numGeo) = TRIM(tag) IF (numGeo < mxGeo) GO TO 10 WRITE (iUnitT, 90) mxGeo 90 FORMAT (/' SUBPROGRAM GetGEO READ ONLY', I6, ' BENCHMARKS.' & & /' INCREASE PARAMETER maxGeo AND RECOMPILE.') CALL Pause() STOP 100 RETURN END SUBROUTINE GetGEO SUBROUTINE GetMOR (iUnitM, iUnitT, mxMOR, & ! INTENT(IN) & MORTag, MORPhi, MORThe, & ! INTENT(OUT) & MORVel, MORSig, numMOR) ! Reads total seafloor spreading rates at mid-ocean rises. ! Current FORMAT matches file "SPREADIN.NUVEL1". ! MORTag = ridge-site name ! MORPhi = East longitude (radians) ! MORThe = co-latitude (radians) ! MORVel = velocity magnitude, m/s ! MORSig = standard deviation of velocity ! IF FIRST READ FAILS (NO FILE CONNECTED), RETURNS numMOR=0 USE DSphere ! in Peter Bird's file DSphere.f90 IMPLICIT NONE INTEGER, INTENT(IN) :: iUnitM, iUnitT, mxMOR CHARACTER*5, DIMENSION(mxMOR), INTENT(OUT) :: MORTag REAL*8, DIMENSION(mxMOR), INTENT(OUT) :: MORPhi, MORThe, MORVel, MORSig INTEGER, INTENT(OUT) :: numMOR CHARACTER*5 :: tag INTEGER :: ios REAL*8 :: eLon, nLat, v, sigma REAL*8 :: secPYr secPYr = 365.25D0 * 24.0D0 * 60.0D0 * 60.0D0 WRITE(iUnitT, 1) iUnitM 1 FORMAT(/' ATTEMPTING TO READ SEAFLOOR SPREADING RATES FROM UNIT', I3/) numMOR = 0 READ (iUnitM, * , IOSTAT = ios) IF (ios /= 0) RETURN READ (iUnitM, * ) READ (iUnitM, * ) 10 READ (iUnitM, 11, END = 100) tag, nLat, eLon, v, sigma 11 FORMAT (A5, F7.2, F8.2, F6.1, F5.1) numMOR = numMOR + 1 MORTag(numMOR) = tag MORPhi(numMOR) = eLon * radians_per_degree MORThe(numMOR) = (90.0D0 - nLat) * radians_per_degree MORVel(numMOR) = v / (1000.0D0 * secPYr) MORSig(numMOR) = sigma / (secPYr * 1000.0D0) IF (numMOR < mxMOR) GO TO 10 WRITE (iUnitT, 90) mxMOR 90 FORMAT (/' SUBPROGRAM GetMOR READ ONLY', I6, ' SPREADING RATES.' & & /' INCREASE PARAMETER maxMOR AND RECOMPILE.') CALL Pause() STOP 100 RETURN END SUBROUTINE GetMOR SUBROUTINE GetNet (iUnit7, iUnit8, & ! input & mxDOF, mxEl, mxFEl, mxNode, & & brief, continuum_LRi, cooling_curvature, & ! output & density_anomaly, & & dQdTdA, elev, fault_LRi, fDip, & & nFakeN, nFl, nodeF, nodes, nRealN, & & numEl, numNod, n1000, offMax, offset, & & title1, tLNode, xNode, yNode, zMNode, & & checkE, checkF, checkN) ! work ! Read finite element grid from unit iUnit7 (assumed already OPENed). ! Echoes the important values to unit iUnit8 (assumed already OPENed). IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - INTEGER, INTENT(IN) :: iUnit7, iUnit8, mxDOF, mxEl, mxFEl, mxNode ! input LOGICAL, INTENT(OUT) :: brief ! output INTEGER, INTENT(OUT) :: continuum_LRi ! output REAL*8, INTENT(OUT) :: cooling_curvature, density_anomaly, dQdTdA, elev ! output INTEGER, INTENT(OUT) :: fault_LRi ! output REAL*8, INTENT(OUT) :: fDip ! output INTEGER, INTENT(OUT) :: nFakeN, nFl, nodeF, nodes, nRealN, numEl, numNod, n1000 ! output REAL*8, INTENT(OUT) :: offMax, offset ! output CHARACTER*80, INTENT(OUT) :: title1 ! output REAL*8, INTENT(OUT) :: tLNode, xNode, yNode, zMNode ! output LOGICAL checkE, checkF, checkN ! work ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - CHARACTER*80 :: longer_line, shorter_line LOGICAL allOK INTEGER i, index, j, k, l, LRi, n, nrt2 REAL*8 cooling_curvature_cpm2, density_anomaly_kgpm3, dips, elevi, & & off, pLat, pLon, qi, tli, vector, xi, yi, zmi DIMENSION checkE(mxEl), checkF(mxFEl), checkN(mxNode), & & continuum_LRi(mxEl), & & cooling_curvature(mxNode), & & density_anomaly(mxNode), & & dQdTdA(mxNode), elev(mxNode), & & fault_LRi(mxFEl), & & fDip(2, mxFEl), nodeF(4, mxFEl), & & nodes(3, mxEl), offset(mxFEl), tLNode(mxNode), & & xNode(mxNode), yNode(mxNode), zMNode(mxNode) DIMENSION dips(3), vector(9) title1 = ' ' READ (iUnit7, 2) title1 2 FORMAT (A80) WRITE (iUnit8, 3) title1 3 FORMAT(/' Title of finite-element grid ='/' ',A80) ! Read number of nodes, plus out-dated parameters that once ! permitted boundary nodes to be specially numbered as ! "fake" nodes with numbers from n1000+1 ... n1000+nFakeN. ! This option is no longer supported by my programs! ! (Option "brief" suppresses most output.) READ (iUnit7, * ) numNod, nRealN, nFakeN, n1000, brief IF (numNod /= (nRealN + nFakeN)) THEN WRITE (iUnit8, 4) numNod, nRealN, nFakeN 4 FORMAT (/' ERR0R: numNod (',I6,') IS NOT EQUAL TO SUM' & & /' OF nRealN (',I6,') AND nFakeN (',I6,').') CALL Pause() STOP END IF IF (nRealN > n1000) THEN WRITE (iUnit8, 5) nRealN, n1000 5 FORMAT (/' ERR0R: nRealN (',I6,') IS GREATER THAN' & & /' n1000 (',I6,').') CALL Pause() STOP END IF IF (numNod > mxNode) THEN WRITE (iUnit8, 10) numNod 10 FORMAT(/' INCREASE ARRAY-SIZE maxNod TO BE AT LEAST' & & /' THE NUMBER OF NODES (', I6, ') AND RECOMPILE.') CALL Pause() STOP END IF nrt2 = nRealN * 2 IF (nrt2 > mxDOF) THEN WRITE (iUnit8, 12) nrt2 12 FORMAT (/' INCREASE ARRAY-SIZE maxDOF TO ', I6, & & ' AND RECOMPILE.') CALL Pause() STOP END IF IF (brief) THEN WRITE (iUnit8, 35) 35 FORMAT(/' (Since option ""brief"" = .TRUE., grid will not be echoed here.)') ELSE WRITE (iUnit8, 40) numNod 40 FORMAT (/' There are',I5,' nodes in the grid') WRITE (iUnit8, 50) 50 FORMAT (/ & & 77X,' mantle'/ & & 77X,' crustal lithosphere'/ & & ' node E-longitude N-latitude', & & ' theta phi elevation', & & ' heat-flow thickness thickness'/) END IF DO 90 k = 1, numNod checkN(k) = .FALSE. 90 CONTINUE DO 100 k = 1, numNod CALL ReadN (iUnit7, iUnit8, 9, & ! input & vector) ! output index = vector(1) + 0.5D0 IF (index > nRealN) THEN IF ((index <= n1000).OR. & & (index > (n1000 + nFakeN))) THEN WRITE (iUnit8, 91) index 91 FORMAT (/' ERR0R: ILLEGAL NODE NUMBER: ',I6) CALL Pause() STOP END IF END IF pLon = vector(2) pLat = vector(3) IF (ABS(pLat) > 90.01) THEN WRITE (iUnit8, 92) index 92 FORMAT (/' ERR0R: ABS(latitude) > 90 AT NODE ',I6) CALL Pause() STOP END IF IF (ABS(pLat) > 89.99D0) THEN WRITE (iUnit8, 93) index 93 FORMAT (/' ERR0R: NODE ',I6,' LIES ON A POLE.' & & /' THIS IS A SINGULAR POINT OF THE' & & ,' SPHERICAL COORDINATE SYSTEM.' & & /' MOVE THIS NODE, AT LEAST SLIGHTLY.') CALL Pause() STOP END IF xi = (90.0D0 - pLat) * 0.0174532925199433D0 yi = pLon * 0.0174532925199433D0 elevi = vector(4) qi = vector(5) zmi = vector(6) tli = vector(7) density_anomaly_kgpm3 = vector(8) cooling_curvature_cpm2 = vector(9) IF (index <= nRealN) THEN i = index ELSE i = nRealN + index - n1000 END IF checkN(i) = .TRUE. xNode(i) = xi yNode(i) = yi elev(i) = elevi dQdTdA(i) = qi IF (qi < 0.0D0) THEN WRITE (iUnit8, 96) 96 FORMAT (' NEGATIVE HEAT-FLOW IS NON-PHYSICAL.') CALL Pause() STOP END IF IF (zmi < 0.0D0) THEN WRITE (iUnit8, 97) 97 FORMAT(' NEGATIVE CRUSTAL THICKNESS IS NON-PHYSICAL.') CALL Pause() STOP END IF zMNode(i) = zmi IF (tli < 0.0D0) THEN WRITE (iUnit8, 98) 98 FORMAT(' NEGATIVE MANTLE LITHOSPHERE THICKNESS IS', & & ' NON-PHYSICAL.') CALL Pause() STOP END IF tLNode(i) = tli IF (.NOT.brief) THEN WRITE (iUnit8, 99) INDEX, pLon, pLat, xi, yi, elevi, & & qi, zmi, tli 99 FORMAT (' ',I10,0P,2F12.3,2F11.5,1P,3E10.2,E12.2) END IF density_anomaly(i) = density_anomaly_kgpm3 cooling_curvature(i) = cooling_curvature_cpm2 100 CONTINUE allOK = .TRUE. DO 101 i = 1, numNod allOK = allOK.AND.checkN(i) 101 CONTINUE IF (.NOT.allOK) THEN WRITE (iUnit8, 102) 102 FORMAT(' THE FOLLOWING NODES WERE NEVER READ:') DO 104 i = 1, numNod IF (i <= nRealN) THEN index = i ELSE index = n1000 + i - nRealN END IF IF (.NOT.checkN(i)) WRITE(iUnit8, 103) index 103 FORMAT (' ',36X,I6) 104 CONTINUE CALL Pause() STOP END IF ! Read triangular elements: READ (iUnit7, * ) numEl IF (numEl > mxEl) THEN WRITE (iUnit8, 108) numEl 108 FORMAT(/' INCREASE PARAMETER maxEl TO BE AT LEAST EQUAL' & & /' TO THE NUMBER OF ELEMENTS (',I6,') AND RECOMPILE.') CALL Pause() STOP END IF DO 109 k = 1, numEl checkE(k) = .FALSE. 109 CONTINUE IF (.NOT.brief) THEN WRITE (iUnit8, 110) numEl 110 FORMAT(/' There are ',I6,' triangular continuum elements.'/ & & ' (Node numbers for each are given at corners, counter', & & 'clockwise'/ / & & ' element C1 C2 C3 Lithospheric_Rheology#') END IF DO 200 k = 1, numEl ! (Elements need not be input in order, but must all be present.) READ (iUnit7, "(A)") longer_line CALL Extract_LRi(longer_line, LRi, shorter_line) continuum_LRi(k) = LRi READ (shorter_line, * ) i, (nodes(j, i), j = 1, 3) IF ((i < 1).OR.(i > numEl)) THEN WRITE (iUnit8, 115) i 115 FORMAT (/' ERR0R: ILLEGAL ELEMENT NUMBER: ',I6) CALL Pause() STOP END IF checkE(i) = .TRUE. IF (.NOT.brief) THEN IF (LRi == 0) THEN WRITE (iUnit8, 120) i, (nodes(j, i), j = 1, 3) 120 FORMAT (' ', I6, ':', 3I10) ELSE WRITE (iUnit8, 121) i, (nodes(j, i), j = 1, 3), LRi 121 FORMAT (' ', I6, ':', 3I10, ' LR', I8) END IF END IF DO 130 j = 1, 3 n = nodes(j, i) IF (n > nRealN) n = nRealN + (n - n1000) IF ((n <= 0).OR.(n > numNod)) THEN WRITE (iUnit8, 125) nodes(j, i) 125 FORMAT (' NODE NUMBER ',I6,' IS ILLEGAL.') CALL Pause() STOP END IF nodes(j, i) = n 130 CONTINUE 200 CONTINUE allOK = .TRUE. DO 201 i = 1, numEl allOK = allOK.AND.checkE(i) 201 CONTINUE IF (.NOT.allOK) THEN WRITE (iUnit8, 202) 202 FORMAT (' THE FOLLOWING ELEMENTS WERE NEVER READ:') DO 204 i = 1, numEl IF (.NOT.checkE(i)) WRITE(iUnit8, 203) i 203 FORMAT (' ',39X,I6) 204 CONTINUE CALL Pause() STOP END IF ! Read fault elements: READ (iUnit7, * ) nFl IF (nFl > mxFEl) THEN WRITE (iUnit8, 220)nFl 220 FORMAT (/' INCREASE PARAMETER maxFEl TO BE AT LEAST EQUAL' & & /' TO THE NUMBER OF FAULTS (',I6,') AND RECOMPILE.') CALL Pause() STOP END IF offMax = 0.0D0 DO 222 i = 1, nFl checkF(i) = .FALSE. 222 CONTINUE IF (.NOT.brief) WRITE(iUnit8, 230) nFl 230 FORMAT(/ /' There are ', I6, ' great-circle fault elements.') IF ((.NOT.brief).AND.(nFl > 0)) WRITE(iUnit8, 231) 231 FORMAT (/' (The 4 node numbers defining each element must be', & & ' in a counterclockwise order:'/ & & ' n1, and n2 are in left-to-right seguence on the', & & ' near side,'/ & & ' then n3 is opposite n2, and n4 is opposite n1.'/, & & ' (Fault dips are given at n1, n2,', & & ' in degrees from horizontal;'/ & & ' positive dips are toward n1 and n2, respectively, '/ & & ' while negative dips are toward n4 and n3.)'/ & & ' Offset is the total past slip of the fault.'/ / & & ' Element n1 n2 n3 n4 dip1 dip2', & & ' offset Lithospheric_Rheology#'/) 240 FORMAT (' ', I6, ':', 4I5, 1X, 2F6.1, 1X, F9.0) DO 300 k = 1, nFl off = 0.0D0 READ(iUnit7, "(A)") longer_line CALL Extract_LRi(longer_line, LRi, shorter_line) fault_LRi(k) = LRi READ(shorter_line, * ) i, (nodeF(j, k), j = 1, 4), (dips(l), l = 1, 2), off IF ((i < 1).OR.(i > nFl)) THEN WRITE (iUnit8, 241) i 241 FORMAT (/' ERR0R: ILLEGAL FAULT NUMBER: ', I6) CALL Pause() STOP END IF checkF(i) = .TRUE. IF (.NOT.brief) THEN IF (LRi == 0) THEN WRITE (iUnit8, 240) i, (nodeF(j, i), j = 1, 4), (dips(l), l = 1, 2), off ELSE WRITE (iUnit8, 242) i, (nodeF(j, i), j = 1, 4), (dips(l), l = 1, 2), off, LRi 242 FORMAT (' ', I6, ':', 4I5, 1X, 2F6.1, 1X, F9.0, " LR", I8) END IF END IF DO 250 j = 1, 4 n = nodeF(j, i) IF (n > nRealN) n = nRealN + (n - n1000) IF ((n <= 0).OR.(n > numNod)) THEN WRITE (iUnit8, 243) nodeF(j, i), i 243 FORMAT (/' ERR0R: ILLEGAL NODE NUMBER (', I6, ') IN FAULT ',I6) CALL Pause() STOP END IF nodeF(j, i) = n 250 CONTINUE DO 260 l = 1, 2 IF (ABS(dips(l)) > 90.0D0) THEN WRITE(iUnit8, 252) dips(l) 252 FORMAT(' ILLEGAL DIP OF ',F10.4,'; SHOULD BE IN', & & ' RANGE OF -90. TO +90. DEGREES.'/ & & ' (NOTE: ALL DIPS ARE IN DEGREES FROM THE', & & ' HORIZONAL;'/ & & ' A + PREFIX (OR NONE) INDICATES A DIP', & & ' TOWARD THE n1-n2 SIDE;'/ & & ' A - PREFIX INDICATES A DIP TOWARD', & & ' THE n4-n3 SIDE.)') CALL Pause() STOP END IF IF (dips(l) < 0.0D0) dips(l) = 180.0D0 + dips(l) fDip(l, i) = dips(l) * 0.0174532925199433D0 260 CONTINUE IF (off < 0.0D0) THEN WRITE (iUnit8, 280) off 280 FORMAT (' ILLEGAL FAULT OFFSET OF ',1P,E10.2, & & ' FOR FAULT ELEMENT',I6/ & & ' OFFSETS MAY NOT BE NEGATIVE.') CALL Pause() STOP END IF offset(i) = off offMax = MAX(offMax, off) 300 CONTINUE allOK = .TRUE. DO 301 i = 1, nFl allOK = allOK.AND.checkF(i) 301 CONTINUE IF (.NOT.allOK) THEN WRITE(iUnit8, 302) 302 FORMAT(' THE FOLLOWING FAULTS WERE NEVER READ:') DO 304 i = 1, nFl IF (.NOT.checkF(i)) WRITE(iUnit8, 303) i 303 FORMAT(' ',36X,I6) 304 CONTINUE CALL Pause() STOP ELSE IF (offMax > 0.0D0) THEN WRITE (iUnit8, 400) offMax 400 FORMAT (/' Greatest fault offset read was ',1P,D10.2) ELSE WRITE (iUnit8, 401) 401 FORMAT (/' Since fault offsets are all zero,', & & ' input parameter Byerly will have no effect.') END IF END IF IF (.NOT. brief) WRITE (iUnit8, 999) 999 FORMAT (' --------------------------------------------------', & & '-----------------------------') END SUBROUTINE GetNet SUBROUTINE GetPBx (iUnitF, iUnitT, names, nPBnd, nPlate, & ! input & nDPlat, pLat, pLon) ! output ! Sets up arrays defining the plates in the PB2002 model of: ! Bird [2003; G**3]. ! (The rotation vectors of the plates are contained in DATA ! statements in the main PROGRAM.) ! The digitized boundaries of the plates (continuous closed curves, ! always circling counterclockwise, and redundantly describing ! each plate boundary *twice*, from each side) ! are read here, from an input file such as 'PB2002_plates.dig', ! on Fortran input device iUnitF. ! The convention for identifying the plates is a 2-character symbol. ! See array "names" in the main PROGRAM. !------------------------------------------------------- IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - INTEGER, INTENT(IN) :: iUnitF, iUnitT ! input CHARACTER*2, INTENT(IN) :: names ! input INTEGER, INTENT(IN) :: nPBnd, nPlate ! input INTEGER, INTENT(OUT) :: nDPlat ! output REAL*8, INTENT(OUT) :: pLat, pLon ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - CHARACTER*2 symbol CHARACTER*3 stars INTEGER i, ios, ip, l, nRead DIMENSION names(nPlate), nDPlat(nPlate) DIMENSION pLat(nPlate, nPBnd), pLon(nPlate, nPBnd) !------------------------------------------------------ nRead = 0 100 READ (iUnitF, 101, END = 201, IOSTAT = ios) symbol IF((nRead == 0).AND.(ios /= 0)) THEN WRITE(*, "(' ERR','OR:File not found, or file empty.')") CALL Pause() STOP END IF 101 FORMAT (A2) DO 120 l = 1, nPlate IF(symbol == names(l)) THEN ip = l GO TO 140 END IF 120 CONTINUE WRITE (iUnitT, 121) symbol, iUnitF 121 FORMAT (/' ERR0R: BAD PLATE NAME ', A, ' ON INPUT DEVICE ', I3) CALL Pause() STOP 140 nRead = nRead + 1 IF (nRead > nPlate) THEN WRITE (iUnitT, 141) WRITE (*, 141) 141 FORMAT(/' Increase nPlate and recompile.') CALL Pause() STOP END IF i = 0 142 READ (iUnitF, 145, END = 201) stars 145 FORMAT (A3) IF (stars == '***') THEN nDPlat(ip) = i GO TO 100 END IF BACKSPACE iUnitF i = i + 1 IF (i > nPBnd) THEN WRITE (iUnitT, 146) WRITE (*, 146) 146 FORMAT(/' Increase nPBnd and recompile.') CALL Pause() STOP END IF READ (iUnitF, * ) pLon(ip, i), pLat(ip, i) pLon(ip, i) = pLon(ip, i) * 0.0174532925199433D0 pLat(ip, i) = pLat(ip, i) * 0.0174532925199433D0 GO TO 142 201 IF (nRead < nPlate) THEN WRITE(iUnitT, "(' ERR','OR: Expecting ', I3, ' plates but' & & /' read outlines of only ', I3)") nPlate, nRead WRITE(*, "(' ERR','OR: Expecting ', I3, ' plates but' & & /' read outlines of only ', I3)") nPlate, nRead CALL Pause() STOP END IF RETURN END SUBROUTINE GetPBx SUBROUTINE GetSLF(iUnitD, iUnitT, mxGSR, & ! INTENT(IN) & fName, rLat, delVZ, fltLon, fltLat, & ! INTENT(OUT) & nFData) !READ in fault slip rates, for scoring purposes. IMPLICIT NONE INTEGER, INTENT(IN) :: iUnitD, iUnitT, mxGSR CHARACTER*100, DIMENSION(mxGSR), INTENT(OUT) :: fName REAL*8, DIMENSION(2, mxGSR), INTENT(OUT) :: rLat, delVZ REAL*8, DIMENSION(mxGSR), INTENT(OUT) :: fltLon, fltLat INTEGER, INTENT(OUT) :: nFData CHARACTER*8 :: c8r1, c8r2, c8v1, c8v2 INTEGER :: i, ios LOGICAL :: dummy ! Detech whether the file exists? nFData = 0 WRITE (iUnitT, 101) iUnitD 101 FORMAT (/' ATTEMPTING TO READ FAULT SLIP RATES FROM UNIT', I3/) READ (iUnitD, * , IOSTAT = ios) IF (ios == 0) THEN REWIND(iUnitD) ELSE RETURN END IF DO 106 i = 1, mxGSR fName(nFData + 1) = ' ' READ(iUnitD, "(A)", END = 111) fName(nFData + 1) ! NOTE: The line with four slip-rate numbers (in mm/year) ! is marked with initial letter "F". ! This is to prevent program -Projector- ! from interpreting such lines as (x, y) data. ! When we read such lines here, we just assign the letter "F" ! to LOGICAL DUMMY (and never use it). READ(iUnitD, * , END = 111) dummy, & & rLat(1, nFData + 1), rLat(2, nFData + 1), & & delVZ(1, nFData + 1), delVZ(2, nFData + 1) READ(iUnitD, * , END = 111) fltLon(nFData + 1), fltLat(nFData + 1) READ(iUnitD, * , END = 111) nFData = nFData + 1 106 CONTINUE WRITE (iUnitT, 107) mxGSR 107 FORMAT (/' ERR0R: PARAMETER maxDat = ', I6, ' IS NOT LARGE'/ & & ' ENOUGH FOR THE NUMBER OF GEOLOGIC SLIP-RATE DATA.') CALL Pause() STOP 111 IF (nFData == 0) THEN WRITE (iUnitT, 112) iUnitD 112 FORMAT (/' NO FAULT SLIP RATE DATA COULD BE READ FROM UNIT ', I2) ELSE WRITE (iUnitT, 113) nFData, iUnitD 113 FORMAT (/' ', I4, ' data on fault slip rates were read from unit ', I2, ':'/) WRITE (iUnitT, 114) !header line of data table; adjust length to agree with FORMAT 145 below. 114 FORMAT ( & & ' ', & & ' -in millimeters/year- -in degrees-'/ & & ' ', & & ' ' & & ' Min. Max. Min. Max. East North'/ & & ' Fault ', & & ' ' & & 'Dextral Dextral Throw Throw Longitude Latitude'/ & & ' ---------------------------------------------------------------------------------------------------- ', & & '------- ------- ------- ------- --------- --------') DO 150 i = 1, nFData IF (rLat(1, i) == 0.0D0) THEN c8r1 = ' ?' ELSE WRITE (c8r1, "(F8.2)") rLat(1, i) END IF IF (rLat(2, i) == 0.0D0) THEN c8r2 = ' ?' ELSE WRITE (c8r2, "(F8.2)") rLat(2, i) END IF IF (delVZ(1, i) == 0.0D0) THEN c8v1 = ' ?' ELSE WRITE (c8v1, "(F8.2)") delVZ(1, i) END IF IF (delVZ(2, i) == 0.0D0) THEN c8v2 = ' ?' ELSE WRITE (c8v2, "(F8.2)") delVZ(2, i) END IF WRITE(iUnitT, 145)fName(i)(1:100), & & c8r1, c8r2, & & c8v1, c8v2, & & fltLon(i), fltLat(i) 145 FORMAT(' ', A, 4A8, F10.3, F9.3) 150 CONTINUE WRITE (iUnitT, 170) 170 FORMAT (' NOTE: Throw-rates are always positive in', & & ' cases of convergence (thrusting),' / & & ' and negative in cases of divergence', & & ' (normal faulting).') END IF END SUBROUTINE GetSLF SUBROUTINE GetSKS (iUnitK, iUnitT, mxSKS, & ! INTENT(IN) & numSKS, SKS_tag, SKS_theta, SKS_phi, & ! INTENT(OUT) & SKS_argument, SKS_delay) ! READs data on upper-mantle seismic anisotropy, from SKS splitting. ! Current FORMAT matches file "Fouch_upper-mantle_SKS_splitting-2004.dat": ! If first READ fails (no file connected), returns numSKS = 0. USE DSphere ! in Peter Bird's file DSphere.f90 IMPLICIT NONE INTEGER, INTENT(IN) :: iUnitK ! Fortran device number to read from INTEGER, INTENT(IN) :: iUnitT ! Fortran devide number for error messages INTEGER, INTENT(IN) :: mxSKS ! maximum number of SKS data that can be stored INTEGER, INTENT(OUT) :: numSKS ! number of SKS data actually read CHARACTER*5, DIMENSION(mxSKS), INTENT(OUT) :: SKS_tag ! 5-byte ID code for datum REAL*8, DIMENSION(mxSKS), INTENT(OUT) :: SKS_theta ! colatitude, in radians REAL*8, DIMENSION(mxSKS), INTENT(OUT) :: SKS_phi ! colatitude, in radians REAL*8, DIMENSION(mxSKS), INTENT(OUT) :: SKS_argument ! phi or fast direction, in radians counterclockwise from S REAL*8, DIMENSION(mxSKS), INTENT(OUT) :: SKS_delay ! SKS splitting time, in s CHARACTER*5 :: tag INTEGER :: azimuth, ios LOGICAL :: problem REAL*8 :: dt, latitude, longitude WRITE(iUnitT, 1) iUnitK 1 FORMAT(/' ATTEMPTING TO READ FAST SKS AZIMUTHS AND DELAYS FROM UNIT', I3/) numSKS = 0 ! Do some test READs to see if a file is connected: problem = .FALSE. READ (iUnitK, "(A)", IOSTAT = ios) tag ! title line? problem = problem.OR.(ios /= 0) READ (iUnitK, "(A)", IOSTAT = ios) tag ! header line? problem = problem.OR.(ios /= 0) IF (problem) RETURN ! Now proceed to read actual datum lines: ! Main loop: 10 READ (iUnitK, 32, END = 100) tag, latitude, longitude, azimuth, dt 32 FORMAT (A5, 1X, 2F10.4, I8, F8.2) numSKS = numSKS + 1 SKS_tag(numSKS) = tag SKS_phi(numSKS) = longitude * radians_per_degree SKS_theta(numSKS) = (90.0D0 - latitude) * radians_per_degree SKS_argument(numSKS) = (180.0D0 - azimuth) * radians_per_degree SKS_delay(numSKS) = dt IF (numSKS < mxSKS) GO TO 10 WRITE (iUnitT, 90) mxSKS 90 FORMAT (/' Subprogram GetSTR READ only', I6, ' data points.' & & /' Increase PARAMETER maxSKS and recompile.') CALL Pause() STOP 100 RETURN END SUBROUTINE GetSKS SUBROUTINE GetSTR (iUnitS, iUnitT, mxStr, & ! INTENT(IN) & strTag, strThe, strPhi, & ! INTENT(OUT) & strArg, strQua, strReg, numStr) ! Reads most-compressive horizontal principal stress azimuths. ! Current FORMAT matches file "WSM1092.LRECL116": ! strTag = short 5-byte identification ! strPhi = East longitude (radians) ! strThe = co-latitude (radians) ! strArg = azimuth of sigme_1H, in radians counterclockwise from South. ! strQua = relative quality, converted from letter to REAL*8 weight. ! strReg = 2-letter abbreviation for stress regime. ! If first READ fails (no file connected), returns numStr = 0. USE DSphere ! in Peter Bird's file DSphere.f90 IMPLICIT NONE INTEGER, INTENT(IN) :: iUnitS, iUnitT, mxStr CHARACTER*5, DIMENSION(mxStr), INTENT(OUT) :: strTag REAL*8, DIMENSION(mxStr), INTENT(OUT) :: strThe, strPhi, strArg, strQua CHARACTER*2, DIMENSION(mxStr), INTENT(OUT) :: strReg INTEGER, INTENT(OUT) :: numStr CHARACTER*5 :: tag CHARACTER*2 :: regime CHARACTER*1 :: qualit REAL*8 :: eLon, nLat INTEGER :: azimut, ios WRITE(iUnitT, 1) iUnitS 1 FORMAT(/' ATTEMPTING TO READ STRESS DIRECTIONS FROM UNIT', I3/) numStr = 0 ! Do one test READ to see if a file is connected? READ (iUnitS, 11, IOSTAT = ios) tag, nLat, eLon, azimut, qualit, regime IF (ios == 0) THEN BACKSPACE iUnitS ELSE RETURN END IF ! Main loop: 10 READ (iUnitS, 11, END = 100) tag, nLat, eLon, azimut, qualit, regime 11 FORMAT (A5, 2F9.3, 7X, I3, 17X, A1, 5X, A2) IF ((azimut.LT.-180.0D0).OR.(azimut.GT.360.0D0)) THEN WRITE (*,"(' ERROR: Bad azimuth in stress data: ', I6)") azimut WRITE (*, 12) tag, nLat, eLon, azimut, qualit, regime 12 FORMAT (' ', A5, 2F9.3, 7X, I3, 17X, A1, 5X, A2) STOP END IF numStr = numStr + 1 strTag(numStr) = tag strPhi(numStr) = eLon * radians_per_degree strThe(numStr) = (90.0D0 - nLat) * radians_per_degree strArg(numStr) = (180.0D0 - azimut) * radians_per_degree IF (qualit == 'A') THEN strQua(numStr) = 4.0D0 ELSE IF (qualit == 'B') THEN strQua(numStr) = 3.0D0 ELSE IF (qualit == 'C') THEN strQua(numStr) = 2.0D0 ELSE IF (qualit == 'D') THEN strQua(numStr) = 1.0D0 ELSE strQua(numStr) = 0.0D0 END IF strReg(numStr) = regime IF (numStr < mxStr) GO TO 10 WRITE (iUnitT, 90) mxStr 90 FORMAT (/' SUBPROGRAM GetSTR READ ONLY', I6, ' AZIMUTHS.' & & /' INCREASE PARAMETER maxStr AND RECOMPILE.') CALL Pause() STOP 100 RETURN END SUBROUTINE GetSTR SUBROUTINE Interp (fAtNod, mxEl, mxNode, nodes, numEl, & ! input & fAtIP) ! output ! Interpolate a scalar function known at the nodes (fAtNod) ! to values at the 7 integration points in each triangular ! continuum element. Note that simple linear interpolation in ! a plane-triangle is used. Thus, this routine is NOT suitable ! for interpolating velocity vectors from nodes to integration points. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - REAL*8, INTENT(IN) :: fAtNod ! input INTEGER, INTENT(IN) :: mxEl, mxNode, nodes, numEl ! input REAL*8, INTENT(OUT) :: fAtIP ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - DOUBLE PRECISION points COMMON / S1S2S3 / points DIMENSION points(3, 7) ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - INTEGER i, m DIMENSION fAtNod(mxNode), fAtIP(7, mxEl), nodes(3, mxEl) DO 100 m = 1, 7 DO 90 i = 1, numEl fAtIP(m, i) = points(1, m) * fAtNod(nodes(1, i)) + & & points(2, m) * fAtNod(nodes(2, i)) + & & points(3, m) * fAtNod(nodes(3, i)) 90 CONTINUE 100 CONTINUE RETURN END SUBROUTINE Interp SUBROUTINE Likely (vMin, vMax, v, & ! INTENT(IN) & prob) ! INTENT(OUT) ! Compute the probability that velocity v is consistent with limits ! vMin and vMax (both of which contain random error). IMPLICIT NONE REAL*8, INTENT(IN) :: vMin, vMax, v REAL*8, INTENT(OUT) :: prob INTEGER :: i LOGICAL :: boxcar, post, point REAL*8 :: b, badata, pg, pgf, pv, sdmodv, sdage, spread, tpim1h, v1, v2, vcent, w, & & d1, d2, db, dw, temp, v1u, v1l, v2u, v2l, vu, vl DATA badata / 0.20 / , sdmodv / 0.10 / , sdage / 0.15 / , tpim1h / 0.398942 / prob = 1.0D0 IF ((vMin == 0.0D0) .AND. (vMax == 0.0D0)) RETURN v1 = vMin IF (vMin == 0.0D0) v1 = -1.0D38 v2 = vMax IF (vMax == 0.0D0) v2 = + 1.0D38 vcent = (v1 + v2) / 2.0D0 spread = ABS((v2 - v1) / (v1 + v2)) boxcar = (spread >= 1.3D0 * sdage).OR.(v1 == -1.0D38).OR.(v2 == 1.0D38) point = spread == 0.0D0 post = (.NOT.boxcar).AND.(.NOT.point) vu = v * (1.0D0 + 2.0D0 * sdmodv) vl = v * (1.0D0 - 2.0D0 * sdmodv) IF(v >= 0.0D0) GO TO 4 temp = vu vu = vl vl = temp 4 v1u = v1 * (1.0D0 + 2.0D0 * sdage) v1l = v1 * (1.0D0 - 2.0D0 * sdage) IF(v1 >= 0.0D0) GO TO 6 temp = v1u v1u = v1l v1l = temp 6 v2u = v2 * (1.0D0 + 2.0D0 * sdage) v2l = v2 * (1.0D0 - 2.0D0 * sdage) IF(v2 >= 0.0D0) GO TO 8 temp = v2u v2u = v2l v2l = temp 8 IF(vu > v1l) GO TO 10 prob = badata RETURN 10 IF(vl < v2u) GO TO 20 prob = badata RETURN 20 IF((vl < v1u).OR.(vu > v2l)) GO TO 30 prob = 1.0D0 RETURN 30 prob = 0.0D0 b = -3.1D0 db = 0.10D0 d1 = ABS(v1) * sdage * 1.41421D0 d2 = ABS(v2) * sdage * 1.41421D0 dw = db * sdmodv * ABS(v) w = v - 31.0D0 * dw DO 40 i = 1, 61 b = b + db pv = tpim1h * EXP(-0.5D0 * b**2) w = w + dw IF (point) pgf = EXP(-MIN(99.9D0, 0.5D0 * ((w - v1) / (sdage * v1))**2)) IF (post) pgf = EXP(-MIN(99.9D0, (0.832D0 * (w - vcent) / (spread * vcent))**2)) IF (boxcar) pgf = MIN(0.5D0 * (1.0D0 + ERF((w - v1) / d1)), 1.0D0, 0.5D0 * (1.0D0 + ERF((v2 - w) / d2))) pgf = MAX(pgf, 1.0D-50) pg = badata + (1. - badata) * pgf 40 prob = prob + pv * pg * db END SUBROUTINE Likely SUBROUTINE Limits (area, detJ, iUnitT, mxEl, numEl, & ! input & okDelV, radius, refStr, sphere, tLInt, & & trHMax, zMoho, & & constr, etaMax, fMuMax, visMax) ! output ! Compute area, mean thickness, and other dimensional parameters ! of the plate, then determine values of stiffness limits needed ! to keep velocity errors down to order okDelV at shear stress ! level refStr. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - REAL*8, INTENT(IN) :: area, detJ ! input INTEGER, INTENT(IN) :: iUnitT, mxEl, numEl ! input REAL*8, INTENT(IN) :: okDelV, radius, refStr ! input LOGICAL, INTENT(IN) :: sphere ! input REAL*8, INTENT(IN) :: tLInt, trHMax, zMoho ! input REAL*8, INTENT(OUT) :: constr, etaMax, fMuMax, visMax ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - DOUBLE PRECISION weight COMMON / WgtVec / weight DIMENSION weight(7) ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - INTEGER i, m, nfault REAL*8 dA, side, thick, totalA, totalV, whole DIMENSION area(mxEl), detJ(7, mxEl), tLInt(7, mxEl), zMoho(7, mxEl) totalA = 0.0D0 totalV = 0.0D0 DO 20 m = 1, 7 DO 10 i = 1, numEl da = area(i) * detJ(m, i) * weight(m) totalA = totalA + da totalV = totalV + da * (zMoho(m, i) + tLInt(m, i)) 10 CONTINUE 20 CONTINUE whole = 4.0D0 * 3.14159265358979D0 * radius**2 IF (totalA > (1.02D0 * whole)) THEN WRITE (iUnitT, 21) totalA, whole 21 FORMAT (/' AREA OF GRID (', ES12.4, ') EXCEEDS' & & /' AREA OF PLANET (', ES12.4, '), WHICH MAKES' & & ,' NO SENSE.' & & /' CHECK GRID FOR ABS(LATITUDE) > 90.' & & /' AND ALSO FOR OVERLAPPING ELEMENTS.') CALL Pause() STOP END IF thick = totalV / totalA IF (sphere) THEN side = radius nfault = 1 ELSE side = SQRT(totalA) nfault = 4 END IF constr = nfault * refStr * thick / okDelV etaMax = refStr * thick / (side * okDelV) etaMax = MIN(etaMax, trHMax / okDelV) fMuMax = nfault * refStr / okDelV visMax = 0.25D0 * refStr * side / okDelV WRITE (iUnitT, 50) totalA, totalV, thick, side, constr, etaMax, & & fMuMax, visMax 50 FORMAT (/ /' Subprogram -Limits- performed dimensional analysis'/ & & ' and estimated necessary stiffness limits to balance'/1P, & & ' the conflicting objectives of accuracy and precision:'/ / & & ' area of model = ',D10.3,' length**2'/ & & ' volume of model = ',D10.3,' length**3'/ & & ' typical thickness = ',D10.3,' length'/ & & ' typical width = ',D10.3,' length'/ & & ' constr (constraint weight) = ',D10.3,' force s/length**2'/ & & ' etaMax (max. basal coupling) = ',D10.3,' force s/length**3'/ & & ' fMuMax (max. fault stiffness) = ',D10.3,' force s/length**3'/ & & ' visMax (max. block viscosity) = ',D10.3,' force s/length**2') END SUBROUTINE Limits SUBROUTINE OLDMohr (aCreep, alphaT, bCreep, Biot, Byerly, & ! input & cCreep, cFric, conduc, constr, dCreep, dQdTdA, & & eCreep, elev, fDip, fFric, fMuMax, & & fPFlt, fArg, gMean, & & mxFEl, mxNode, nFl, nodeF, & & offMax, offset, & & oneKm, radio, rhoH2O, rhoBar, & & slide, tauMax, & & tLNode, tSurf, v, wedge, & & zMNode, & & zTranF, & ! modify & fC, fIMuDZ, fPeakS, fSlips, fTStar) ! output ! This subprogram contains the nonlinear rheology of the faults. ! For each of 7 integration points along the length of each fault ! element, it: ! (1) Computes the slip-rate vector on the fault surface; ! (2) Determines the shear stress on the fault surface by Mohr/ ! Coulomb/Navier theory; (This stress is proportional to depth, ! so the calculation is actually done at unit depth and then ! scaled.) ! (3) Proceeds down the dip of the fault, checking temperature, ! strain rate, and pressure to see if frictional or creep ! shear stress is lower; ! (4) Reports the vertical integral of "mu" (the ratio of shear ! stress to slip rate) down the fault as fIMuDZ; ! (Note that the integral is vertical, not on a slant, even though ! conditions are evaluated along a slant path.) ! (5) For dipping, oblique-slip faults only, also reports recommended ! tactical values for the matrix fC and the vector fTStar ! which jointly describe a linearized rheology stiffer than ! the actual nonlinear rheology; ! (6) zTranF is the latest estimate of the depths (in crust, in mantle lithosphere) ! to the brittle/ductile transitions, at the fault midpoint; ! (7) LOGICAL variable fSlips indicates whether the fault is ! slipping at its midpoint; (Otherwise, it is in the artificial ! linearized regime, with stiffness fMuMax.) ! (8) fPeakS gives the peak shear stress at the midpoint of each ! fault, evaluated at the brittle/ductile transition; ! (9) Faults with dip less than "slide" (radians) are limited ! to a maximum down-dip integral shear traction of tauMax. ! Note that pore pressures are considered in the calculation of ! frictional strength: ! *Normal pore pressures reduce the effective normal stress on the ! fault surface by the amount ! -Biot * gMean * rhoH20 * z ! *IF (offMax > 0.) THEN the remaining effective frictional strength ! of the fault is multiplied by the reducing factor ! *(1. - Byerly * offset(i) / offMax). ! This may also be a pore pressure effect, because Byerlee's model is ! that gouge layers have thickness in proportion to offset, and ! that they support non-Darcy static pore pressure gradients which ! allow elevated pore pressures in the core of the gouge, which ! reduce the effective friction of the fault. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - REAL*8, INTENT(IN) :: aCreep, alphaT, bCreep, Biot, Byerly, cCreep, cFric, conduc, & ! input & constr, dCreep, dQdTdA, eCreep, elev, fDip, fFric, fMuMax, & ! input & fPFlt, fArg, gMean ! input INTEGER, INTENT(IN) :: mxFEl, mxNode, nFl, nodeF ! input REAL*8, INTENT(IN) :: offMax, offset, oneKm, radio, rhoH2O, rhoBar, slide ! input REAL*8, INTENT(IN) :: tauMax, tLNode, tSurf ! input DOUBLE PRECISION, INTENT(IN) :: v ! input REAL*8, INTENT(IN) :: wedge, zMNode ! input REAL*8, INTENT(INOUT) :: zTranF ! modify REAL*8, INTENT(OUT) :: fC, fIMuDZ, fPeakS ! output LOGICAL, INTENT(OUT) :: fSlips ! output REAL*8, INTENT(OUT) :: fTStar ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - DOUBLE PRECISION fPhi COMMON / FPhis / fPhi DIMENSION fPhi(4, 7) ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ! Following PARAMETER gives number of steps in vertical integral ! of creep shear stress on ductile parts of faults: INTEGER, PARAMETER :: nStep = 30 ! Higher values obviously cost more. On the other hand, small values ! do not merely approximate the creep law; they also introduce ! some random error which can put a floor on convergence ! of the whole global velocity field. INTEGER i, j, kIter, layer, limit, m, n1, n2, n3, n4 REAL*8 angle, azfull, azhalf, baseZ, cGamma, close, cosr, crust, & & dDPNdZ, delVx, delVy, delTau, dEPdST, dip, dippy, dLEPdC, dLEPdZ, dPMax, dSFdZ, dz, & & efull, ehalf, elevat, ePMoho, fric, huge, mantle, & & normal, oldsc, q, q0, rake, rho, rhoC, & & scfull, schalf, sf0, sFMoho, shearc, shearf, shearp, sheart, sinist, sinr, slip, spread, strain, sum, & & t, t0, t90pc, tand, tAsth, tfull, thalf, thick, thrust, tiny, tlan, tMean, tMeanC, tMoho, topZ, ts, tTrans, tu, tune, & & unitAx, unitAy, unitBx, unitBy, & & vIMuDZ, vitdz, vUpDip, whalf, wfull, z, z0, zAbs, zfull, zhalf, zman, zp, zTrans DOUBLE PRECISION dpt1 LOGICAL locked, pureSS, sloped ! DIMENSIONs of internal convenience arrays: DIMENSION dLEPdZ(2), dSFdZ(2), rho(2), sheart(2), tMean(2), zTrans(2) ! DIMENSIONs of external argument arrays: DIMENSION aCreep(2), alphaT(2), bCreep(2), cCreep(2), conduc(2), & & dCreep(2), dQdTdA(mxNode), elev(mxNode), & & fC(2, 2, 7, mxFEl), fDip(2, mxFEl), & & fIMuDZ(7, mxFEl), fPeakS(2, mxFEl), & & fPFlt(2, 2, 2, 7, mxFEl), fSlips(mxFEl), & & fArg(2, mxFEl), fTStar(2, 7, mxFEl), nodeF(4, mxFEl), & & offset(mxFEl), radio(2), rhoBar(2), & & tauMax(2), tLNode(mxNode), & & v(2, mxNode), zMNode(mxNode), zTranF(2, mxFEl) ! Following two numbers are "very small" and "very large", but not ! so extreme as to cause underflow or overflow. They may need to ! be adjusted, depending on the computer and compiler you use: DATA tiny / 2.0D-38 / DATA huge / 1.0D+38 / cgamma = (1.0D0 + SIN(ATAN(cFric))) / (1.0D0 - SIN(ATAN(cFric))) DO 100 i = 1, nFl IF (offMax <= 0.) THEN fric = fFric ELSE fric = fFric * (1.0D0 - Byerly * offset(i) / offMax) END IF n1 = nodeF(1, i) n2 = nodeF(2, i) n3 = nodeF(3, i) n4 = nodeF(4, i) ! Is this a purely strike-slip fault element? pureSS = (ABS(fDip(1, i) - 1.57079632679490D0) <= wedge).AND. & & (ABS(fDip(2, i) - 1.57079632679490D0) <= wedge) ! If so, compute estimate of relative normal stress ! (relative to vertical stress) by using amount of divergence ! between average of node n1 and n2 and average of node n3 ! and node n4 (in spite on constraint equation): IF (pureSS) THEN !CCCC angle = 0.5 * (fArg(1, i) + fArg(2, i)) !CCCC Line above was replaced due to cycle-shift problem! angle = Chord(fArg(1, i), 0.5D0, fArg(2, i)) unitBx = SIN(angle) unitBy = -COS(angle) ! unitB points outward on the n1-n2 side (away from ! the n3-n4 side). delVx = v(1, n1) * fPFlt(1, 1, 1, 4, i) + v(2, n1) * fPFlt(2, 1, 1, 4, i) & & + v(1, n2) * fPFlt(1, 1, 2, 4, i) + v(2, n2) * fPFlt(2, 1, 2, 4, i) & & - v(1, n3) * fPFlt(1, 1, 2, 4, i) - v(2, n3) * fPFlt(2, 1, 2, 4, i) & & - v(1, n4) * fPFlt(1, 1, 1, 4, i) - v(2, n4) * fPFlt(2, 1, 1, 4, i) delVy = v(1, n1) * fPFlt(1, 2, 1, 4, i) + v(2, n1) * fPFlt(2, 2, 1, 4, i) & & + v(1, n2) * fPFlt(1, 2, 2, 4, i) + v(2, n2) * fPFlt(2, 2, 2, 4, i) & & - v(1, n3) * fPFlt(1, 2, 2, 4, i) - v(2, n3) * fPFlt(2, 2, 2, 4, i) & & - v(1, n4) * fPFlt(1, 2, 1, 4, i) - v(2, n4) * fPFlt(2, 2, 1, 4, i) ! delVx and delVy are the velocities of the n1-n2 side ! relative to the n3-n4 side. spread = delVx * unitBx + delVy * unitBy delTau = constr * spread tlan = 0.5D0 * (tLNode(n1) + tLNode(n2)) zman = 0.5D0 * (zMNode(n1) + zMNode(n2)) IF ((tlan <= 0.0D0).OR.(zTranF(2, i) <= 0.0D0)) THEN ! Crust alone resists convergence: dPMax = -2.0D0 * deltau / zTranF(1, i) dDPNdZ = dpmax / zTranF(1, i) ELSE ! Mantle lithosphere helps to resist convergence: dDPNdZ = -deltau / & & (0.50D0 * zTranF(1, i)**2 + zTranF(2, i) * zman + & & 0.50D0 * zTranF(2, i)**2) END IF ! dDPNdZ is the gradient of excess normal pressure (in ! excess of vertical pressure) with depth on this fault; ! check that it lies within frictional limits of blocks: q = 0.250D0 * (dQdTdA(n1) + dQdTdA(n2) + & & dQdTdA(n3) + dQdTdA(n4)) tTrans = tSurf + zTranF(1, i) * q / conduc(1) - & & zTranF(1, i)**2 * radio(1) / (2. * conduc(1)) tMeanC = (tSurf + tTrans) / 2.0D0 rhoC = rhoBar(1) * (1.0D0 - alphaT(1) * tMeanC) dLEPdC = gMean * (rhoC - rhoH2O * Biot) thrust = dLEPdC * cGamma normal = dLEPdC / cGamma dDPNdZ = MAX(dDPNdZ, normal - dLEPdC) dDPNdZ = MIN(dDPNdZ, thrust - dLEPdC) ELSE ! Different logic will be used; this parameter is not ! really needed. Zero it just to be careful. dDPNdZ = 0.0D0 END IF DO 90 m = 1, 7 ! elevation: elevat = elev(n1) * fPhi(1, m) + elev(n2) * fPhi(2, m) ! heat flow: q = dQdTdA(n1) * fPhi(1, m) + dQdTdA(n2) * fPhi(2, m) ! crustal thickness: crust = zMNode(n1) * fPhi(1, m) + zMNode(n2) * fPhi(2, m) ! mantle lithosphere thickness: mantle = tLNode(n1) * fPhi(1, m) + tLNode(n2) * fPhi(2, m) mantle = MAX(mantle, 0.0D0) ! Moho temperature: tMoho = tSurf + crust * q / conduc(1) - & & crust**2 * radio(1) / (2.0D0 * conduc(1)) ! Temperature at base of plate: tAsth = tMoho + mantle * (q - crust * radio(1)) / conduc(2) - & & mantle**2 * radio(2) / (2.0D0 * conduc(2)) ! mean temperatures: tMean(1) = (tSurf + tMoho) / 2.0D0 tMean(2) = (tMoho + tAsth) / 2.0D0 ! mean densities: rho(1) = rhoBar(1) * (1.0D0 - alphaT(1) * tMean(1)) rho(2) = rhoBar(2) * (1.0D0 - alphaT(2) * tMean(2)) ! derivitives of lithostatic effective pressure wrt depth dLEPdZ(1) = gMean * (rho(1) - rhoH2O * Biot) ePMoho = dLEPdZ(1) * crust dLEPdZ(2) = gMean * (rho(2) - rhoH2O * Biot) ! "angle" is the fault strike, in radians cclkws from +X. !CCCC angle = fArg(1, i) * fPhi(1, m) + fArg(2, i) * fPhi(2, m) !CCCC Line above was replaced due to cycle-shift problem! angle = Chord(fArg(1, i), fPhi(2, m), fArg(2, i)) ! unitA is a unit vector along the fault, from n1 to n2. unitAx = COS(angle) unitAy = SIN(angle) ! unitB is a perpendicular unit vector, pointing out ! toward the n4-n3 side. unitBx = -unitAy unitBy = + unitAx ! relative velocities are for n1-n2 side relative to ! the n4-n3 side: delVx = v(1, n1) * fPFlt(1, 1, 1, m, i) + v(2, n1) * fPFlt(2, 1, 1, m, i) & & + v(1, n2) * fPFlt(1, 1, 2, m, i) + v(2, n2) * fPFlt(2, 1, 2, m, i) & & - v(1, n3) * fPFlt(1, 1, 2, m, i) - v(2, n3) * fPFlt(2, 1, 2, m, i) & & - v(1, n4) * fPFlt(1, 1, 1, m, i) - v(2, n4) * fPFlt(2, 1, 1, m, i) delVy = v(1, n1) * fPFlt(1, 2, 1, m, i) + v(2, n1) * fPFlt(2, 2, 1, m, i) & & + v(1, n2) * fPFlt(1, 2, 2, m, i) + v(2, n2) * fPFlt(2, 2, 2, m, i) & & - v(1, n3) * fPFlt(1, 2, 2, m, i) - v(2, n3) * fPFlt(2, 2, 2, m, i) & & - v(1, n4) * fPFlt(1, 2, 1, m, i) - v(2, n4) * fPFlt(2, 2, 1, m, i) ! sinistral strike-slip rate component: sinist = delVx * unitAx + delVy * unitAy ! convergence rate component (in horizontal plane): close = delVx * unitBx + delVy * unitBy ! dip of the fault (from horizontal on the n1-n2 side): dip = fDip(1, i) * fPhi(1, m) + fDip(2, i) * fPhi(2, m) sloped = ABS(dip - 1.57079632679490D0) > wedge IF (.NOT.sloped) THEN ! case of a near-vertical fault: dSFdZ(1) = (dLEPdZ(1) + dDPNdZ) * fric sFMoho = dSFdZ(1) * crust dSFdZ(2) = (dLEPdZ(2) + dDPNdZ) * fric slip = ABS(sinist) locked = .FALSE. ELSE ! case of a shallow-dipping fault: ! vUpDip is the up-dip velocity component, in the ! fault plane, of the block on the n1-n3 side: vUpDip = close / COS(dip) ! "rake" angle is measured counterclockwise in ! fault plane from horizontal & parallel to "angle": rake = ATan2F(vUpDip, sinist) ! derivitive of effective normal pressure ! with respect to shear traction on fault: dEPdST = TAN(dip) * SIN(rake) ! (Notice that when sense of dip reverses, sign ! change caused by TAN(dip) is cancelled by sign ! change caused by SIN(rake).) ! According to theory, the equation to solve is: ! d(shear traction)/dz = ! fric*(dLEPdZ+dEPdST*d(shear_traction)/dz) ! This may have a physical solution (one with ! positive shear traction). If not, the fault ! is locked. locked = (fric * dEPdST) >= 1.00D0 IF (locked) THEN dSFdZ(1) = huge dSFdZ(2) = huge ELSE dSFdZ(1) = fric * dLEPdZ(1) / (1.00D0 - fric * dEPdST) sFMoho = dSFdZ(1) * crust dSFdZ(2) = fric * dLEPdZ(2) / (1.00D0 - fric * dEPdST) END IF slip = SQRT(sinist**2 + vUpDip**2) END IF slip = MAX(slip, 1.0D8 * tiny) ! Locate plastic/creep transition(s) ! by iterated halving of domain: IF (mantle > 0.) THEN limit = 2 ELSE limit = 1 zTrans(2) = 0.0D0 sheart(2) = 0.0D0 END IF DO 60 layer = 1, limit topZ = 0.0D0 IF (layer == 1) THEN baseZ = crust sf0 = 0.0D0 t0 = tSurf q0 = q z0 = 0.0D0 ELSE baseZ = mantle sf0 = sFMoho t0 = tMoho q0 = q - crust * radio(1) z0 = crust END IF DO 50 kIter = 1, 15 z = 0.50D0 * (topZ + baseZ) zAbs = z + z0 shearf = z * dSFdZ(layer) + sf0 shearp = MIN(shearf, dCreep(layer)) t = t0 + q0 * z / conduc(layer) - (radio(layer) / & & (2.0D0 * conduc(layer))) * z**2 IF (zAbs <= (15.0D0 * oneKm)) THEN t90pc = 0.50D0 * zAbs ELSE IF (zAbs < (45.0D0 * oneKm)) THEN t90pc = (405.0D0 / 8.0D0) * oneKm + & & (-7.0D0) * zAbs + & & (13.0D0 / 40.0D0) * oneKm * (zAbs / oneKm)**2 + & & (-1.0D0 / 300.0D0) * oneKm * (zAbs / oneKm)**3 ELSE t90pc = 2.0D0 * zAbs END IF ! See Turcotte et al (1980) JGR, 85, B11, 6224-6230. strain = slip / t90pc shearc = aCreep(layer) * (strain**eCreep) * & & EXP((bCreep(layer) + cCreep(layer) * z) / t) IF (shearc < shearp) THEN baseZ = z ELSE topZ = z END IF 50 CONTINUE zTrans(layer) = 0.50D0 * (topZ + baseZ) sheart(layer) = zTrans(layer) * dSFdZ(layer) + sf0 60 CONTINUE ! plastic part of vertical integral(s) of traction: ! (A) crust: IF (sheart(1) <= dCreep(1)) THEN vitdz = 0.50D0 * sheart(1) * zTrans(1) ELSE zp = zTrans(1) * dCreep(1) / sheart(1) vitdz = dCreep(1) * (zTrans(1) - 0.50D0 * zp) END IF ! (B) mantle lithosphere: IF ((mantle > 0.).AND.(sheart(2) > sfmoho)) THEN IF (sheart(2) <= dCreep(2)) THEN vitdz = vitdz + 0.50D0 * (sfmoho + sheart(2)) * zTrans(2) ELSE zp = zTrans(2) * (dCreep(2) - sfmoho) / & & (sheart(2) - sfmoho) zp = MAX(zp, 0.) vitdz = vitdz + 0.50D0 * (sfmoho + sheart(2)) * zp + & & dCreep(2) * (zTrans(2) - zp) END IF END IF ! Add creep part(s) of integral, using parabolic rule: sum = 0.0D0 DO 80 layer = 1, limit IF (layer == 1) THEN thick = crust t0 = tSurf q0 = q zAbs = 0.0D0 ELSE thick = mantle t0 = tMoho q0 = q - crust * radio(1) zAbs = crust END IF dz = (thick - zTrans(layer)) / nStep oldsc = sheart(layer) oldsc = MIN(oldsc, dCreep(layer)) z0 = zTrans(layer) DO 70 j = 1, nStep zhalf = z0 + 0.50D0 * dz zfull = z0 + dz azhalf = zhalf + zAbs azfull = zfull + zAbs thalf = t0 + q0 * zhalf / conduc(layer) - & & (radio(layer) / & & (2.0D0 * conduc(layer))) * zhalf**2 tfull = t0 + q0 * zfull / conduc(layer) - & & (radio(layer) / & & (2.0D0 * conduc(layer))) * zfull**2 IF (azhalf <= (15.0D0 * oneKm)) THEN whalf = 0.50D0 * azhalf ELSE IF (azhalf < (45.0D0 * oneKm)) THEN whalf = (405.0D0 / 8.0D0) * oneKm + & & (-7.0D0) * azhalf + & & (13.0D0 / 40.0D0) * oneKm * (azhalf / oneKm)**2 + & & (-1.0D0 / 300.0D0) * oneKm * (azhalf / oneKm)**3 ELSE whalf = 2.0D0 * azhalf END IF IF (azfull <= (15.0D0 * oneKm)) THEN wfull = 0.50D0 * azfull ELSE IF (azfull < (45.0D0 * oneKm)) THEN wfull = (405.0D0 / 8.0D0) * oneKm + & & (-7.0D0) * azfull + & & (13.0D0 / 40.0D0) * oneKm * (azfull / oneKm)**2 + & & (-1.0D0 / 300.0D0) * oneKm * (azfull / oneKm)**3 ELSE wfull = 2.0D0 * azhalf END IF ! See Turcotte et al (1980) JGR, 85, B11, 6224-6230. ehalf = slip / whalf efull = slip / wfull schalf = aCreep(layer) * (ehalf**eCreep) * & & EXP((bCreep(layer) + cCreep(layer) * zhalf) & & / thalf) schalf = MIN(schalf, dCreep(layer)) scfull = aCreep(layer) * (efull**eCreep) * & & EXP((bCreep(layer) + cCreep(layer) * zfull) & & / tfull) scfull = MIN(scfull, dCreep(layer)) sum = sum + dz * (0.1666667D0 * oldsc + & & 0.6666667D0 * schalf + & & 0.1666666D0 * scfull) z0 = zfull oldsc = scfull 70 CONTINUE 80 CONTINUE vitdz = vitdz + sum ! Limit shear traction on subduction zones only: dippy = MIN(dip, 3.14159265358979D0 - dip) IF (dippy <= slide) THEN IF (elevat < 0.0D0) THEN ! apply oceanic subduction zone limit: vitdz = MIN(vitdz, tauMax(1) * SIN(dip)) ELSE ! apply continental subduction zone limit: vitdz = MIN(vitdz, tauMax(2) * SIN(dip)) END IF END IF dpt1 = (1.0D0 * vitdz) / slip vIMuDZ = MIN(dpt1, 1.0D38) fIMuDZ(m, i) = MIN(vimudz, fMuMax * (crust + mantle)) ! Dipping, oblique-slip integration ! points are also characterized ! by fC and fTStar: IF (sloped) THEN ts = sinist * fIMuDZ(m, i) tu = vUpDip * fIMuDZ(m, i) IF (locked) THEN fC(1, 1, m, i) = fIMuDZ(m, i) fC(1, 2, m, i) = 0.0D0 fC(2, 1, m, i) = 0.0D0 fC(2, 2, m, i) = fIMuDZ(m, i) ELSE sinr = SIN(rake) cosr = COS(rake) tand = TAN(dip) ! *** IMPORTANT NOTE: *** ! The following 7 statements are *not* the ! result of theory, but a tactical choice ! which attempts to compromise between ! stability of the linear system, stability ! of the iteration, and efficiency. ! They may be changed if the program does ! no converge satisfactorily! tune = 2.0D0 fC(1, 1, m, i) = fIMuDZ(m, i) * & & (1.0D0 - tune * sinr * cosr**2 * tand) fC(1, 2, m, i) = fIMuDZ(m, i) * & & (tune * cosr**3 * tand) fC(2, 1, m, i) = fIMuDZ(m, i) * & & (-tune * sinr**2 * cosr * tand) fC(2, 2, m, i) = fIMuDZ(m, i) * & & (1.0D0 + tune * sinr * cosr**2 * tand) ! (Often, fC(1,2) is the biggest term. ! In some cases, diagonals become negative. ! For stability, be sure that the fC ! matrix remains positive-definite: fC(1, 1, m, i) = MAX(fC(1, 1, m, i), ABS(fC(1, 2, m, i))) fC(2, 2, m, i) = MAX(fC(2, 2, m, i), ABS(fC(1, 2, m, i))) END IF fTStar(1, m, i) = ts - fC(1, 1, m, i) * sinist - & & fC(1, 2, m, i) * vUpDip fTStar(2, m, i) = tu - fC(2, 1, m, i) * sinist - & & fC(2, 2, m, i) * vUpDip END IF ! Provide interesting diagnostic data at midpoints only: IF (m == 4) THEN fSlips(i) = (.NOT.locked).AND. & & (fIMuDZ(m, i) < (0.99D0 * fMuMax * (crust + mantle))) zTranF(1, i) = zTrans(1) fPeakS(1, i) = MIN(sheart(1), dCreep(1)) zTranF(2, i) = zTrans(2) fPeakS(2, i) = MIN(sheart(2), dCreep(2)) END IF 90 CONTINUE 100 CONTINUE END SUBROUTINE OLDMohr SUBROUTINE Mohr (alphaT, conduc, constr, & ! input & continuum_LRi, & & dQdTdA, elev, & & fault_LRi, fDip, fMuMax, & & fPFlt, fArg, gMean, & & LRn, LR_set_fFric, LR_set_cFric, LR_set_Biot, LR_set_Byerly, & & LR_set_aCreep, LR_set_bCreep, LR_set_cCreep, LR_set_dCreep, LR_set_eCreep, & & mxEl, mxFEl, mxNode, nFl, nodeF, & & offMax, offset, & & oneKm, radio, rhoH2O, rhoBar, & & slide, tauMax, & & tLNode, tSurf, v, wedge, & & zMNode, & & zTranF, & ! modify & fC, fIMuDZ, fPeakS, fSlips, fTStar) ! output ! This subprogram contains the nonlinear rheology of the faults. ! For each of 7 integration points along the length of each fault ! element, it: ! (1) Computes the slip-rate vector on the fault surface; ! (2) Determines the shear stress on the fault surface by Mohr/ ! Coulomb/Navier theory; (This stress is proportional to depth, ! so the calculation is actually done at unit depth and then ! scaled.) ! (3) Proceeds down the dip of the fault, checking temperature, ! strain rate, and pressure to see if frictional or creep ! shear stress is lower; ! (4) Reports the vertical integral of "mu" (the ratio of shear ! stress to slip rate) down the fault as fIMuDZ; ! (Note that the integral is vertical, not on a slant, even though ! conditions are evaluated along a slant path.) ! (5) For dipping, oblique-slip faults only, also reports recommended ! tactical values for the matrix fC and the vector fTStar ! which jointly describe a linearized rheology stiffer than ! the actual nonlinear rheology; ! (6) zTranF is the latest estimate of the depth ! to the brittle/ductile transition, at the fault midpoint; ! (7) LOGICAL variable fSlips indicates whether the fault is ! slipping at its midpoint; (Otherwise, it is in the artificial ! linearized regime, with stiffness fMuMax.) ! (8) fPeakS gives the peak shear stress at the midpoint of each ! fault, evaluated at the brittle/ductile transition; ! (9) Faults with dip less than "slide" (radians) are limited ! to a maximum down-dip integral shear traction of tauMax. ! Note that pore pressures are considered in the calculation of ! frictional strength: ! *Normal pore pressures reduce the effective normal stress on the ! fault surface by the amount ! -Biot * gMean * rhoH20 * z ! *IF (offMax > 0.) THEN the remaining effective frictional strength ! of the fault is multiplied by the reducing factor ! *(1. - Byerly * offset(i) / offMax). ! This may also be a pore pressure effect, because Byerlee's model is ! that gouge layers have thickness in proportion to offset, and ! that they support non-Darcy static pore pressure gradients which ! allow elevated pore pressures in the core of the gouge, which ! reduce the effective friction of the fault. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - REAL*8, INTENT(IN) :: alphaT, conduc, constr ! input INTEGER, INTENT(IN) :: continuum_LRi ! input REAL*8, INTENT(IN) :: dQdTdA, elev ! input INTEGER, INTENT(IN) :: fault_LRi ! input REAL*8, INTENT(IN) :: fDip, fMuMax, fPFlt, fArg, gMean ! input INTEGER, INTENT(IN) :: LRn ! input REAL*8, INTENT(IN) :: LR_set_fFric, LR_set_cFric, LR_set_Biot, LR_set_Byerly, & ! input & LR_set_aCreep, LR_set_bCreep, LR_set_cCreep, LR_set_dCreep, LR_set_eCreep ! input INTEGER, INTENT(IN) :: mxEl, mxFEl, mxNode, nFl, nodeF ! input REAL*8, INTENT(IN) :: offMax, offset, oneKm, radio, rhoH2O, rhoBar, slide ! input REAL*8, INTENT(IN) :: tauMax, tLNode, tSurf ! input DOUBLE PRECISION, INTENT(IN) :: v ! input REAL*8, INTENT(IN) :: wedge, zMNode ! input REAL*8, INTENT(INOUT) :: zTranF ! modify REAL*8, INTENT(OUT) :: fC, fIMuDZ, fPeakS ! output LOGICAL, INTENT(OUT) :: fSlips ! output REAL*8, INTENT(OUT) :: fTStar ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - DOUBLE PRECISION fPhi COMMON / FPhis / fPhi DIMENSION fPhi(4, 7) ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ! Following PARAMETER gives number of steps in vertical integral ! of creep shear stress on ductile parts of faults: INTEGER, PARAMETER :: nStep = 30 ! Higher values obviously cost more. On the other hand, small values ! do not merely approximate the creep law; they also introduce ! some random error which can put a floor on convergence ! of the whole global velocity field. INTEGER i, j, kIter, layer, limit, m, n1, n2, n3, n4 REAL*8 angle, azfull, azhalf, baseZ, cGamma, close, cosr, crust, & & dDPNdZ, delVx, delVy, delTau, dEPdST, dip, dippy, dLEPdC, dLEPdZ, dPMax, dSFdZ, dz, & & efull, ehalf, elevat, ePMoho, fric, huge, mantle, & & normal, oldsc, q, q0, rake, rho, rhoC, & & scfull, schalf, sf0, sFMoho, shearc, shearf, shearp, sheart, sinist, sinr, slip, spread, strain, sum, & & t, t0, t90pc, tand, tAsth, tfull, thalf, thick, thrust, tiny, tlan, tMean, tMeanC, tMoho, topZ, ts, tTrans, tu, tune, & & unitAx, unitAy, unitBx, unitBy, & & vIMuDZ, vitdz, vUpDip, whalf, wfull, z, z0, zAbs, zfull, zhalf, zman, zp, zTrans DOUBLE PRECISION dpt1 LOGICAL locked, pureSS, sloped ! DIMENSIONs of external argument arrays: DIMENSION alphaT(2), conduc(2), & & continuum_LRi(mxEl), & & dQdTdA(mxNode), elev(mxNode), & & fault_LRi(mxFEl), & & fC(2, 2, 7, mxFEl), fDip(2, mxFEl), & & fIMuDZ(7, mxFEl), fPeakS(2, mxFEl), & & fPFlt(2, 2, 2, 7, mxFEl), fSlips(mxFEl), & & fArg(2, mxFEl), fTStar(2, 7, mxFEl), & & LR_set_fFric(0:LRn), LR_set_cFric(0:LRn), LR_set_Biot(0:LRn), LR_set_Byerly(0:LRn), & & LR_set_aCreep(1:2, 0:LRn), LR_set_bCreep(1:2, 0:LRn), LR_set_cCreep(1:2, 0:LRn), LR_set_dCreep(1:2, 0:LRn), LR_set_eCreep(0:LRn), & & nodeF(4, mxFEl), & & offset(mxFEl), radio(2), rhoBar(2), & & tauMax(2), tLNode(mxNode), & & v(2, mxNode), zMNode(mxNode), zTranF(2, mxFEl) ! DECLARATIONS and DIMENSIONs of internal convenience arrays: DIMENSION dLEPdZ(2), dSFdZ(2), rho(2), sheart(2), tMean(2), zTrans(2) INTEGER LRi REAL*8 t_fFric, t_cFric, t_Biot, t_Byerly, t_aCreep(2), t_bCreep(2), t_cCreep(2), t_dCreep(2), t_eCreep ! Following two numbers are "very small" and "very large", but not ! so extreme as to cause underflow or overflow. They may need to ! be adjusted, depending on the computer and compiler you use: DATA tiny / 2.D-38 / DATA huge / 1.D+38 / !Use default rheology to define the environment surrounding all faults: cGamma = (1.0D0 + SIN(ATAN(LR_set_cFric(0)))) / (1.0D0 - SIN(ATAN(LR_set_cFric(0)))) DO 100 i = 1, nFl !Extract the desired rheology for this fault element: LRi = fault_LRi(i) t_fFric = LR_set_fFric(LRi) t_Biot = LR_set_Biot(LRi) t_Byerly = LR_set_Byerly(LRi) t_aCreep(1:2) = LR_set_aCreep(1:2, LRi) t_bCreep(1:2) = LR_set_bCreep(1:2, LRi) t_cCreep(1:2) = LR_set_cCreep(1:2, LRi) t_dCreep(1:2) = LR_set_dCreep(1:2, LRi) t_eCreep = LR_set_eCreep(LRi) !- - - - - - - - - - - - - - - - - - - - - - -- - - - IF (offMax <= 0.0D0) THEN fric = t_fFric ELSE fric = t_fFric * (1.0D0 - t_Byerly * offset(i) / offMax) END IF n1 = nodeF(1, i) n2 = nodeF(2, i) n3 = nodeF(3, i) n4 = nodeF(4, i) ! Is this a purely strike-slip fault element? pureSS = (ABS(fDip(1, i) - 1.57079632679490D0) <= wedge).AND. & & (ABS(fDip(2, i) - 1.57079632679490D0) <= wedge) ! If so, compute estimate of relative normal stress ! (relative to vertical stress) by using amount of divergence ! between average of node n1 and n2 and average of node n3 ! and node n4 (in spite on constraint equation): IF (pureSS) THEN !CCCC angle = 0.5 * (fArg(1, i) + fArg(2, i)) !CCCC Line above was replaced due to cycle-shift problem! angle = Chord(fArg(1, i), 0.5D0, fArg(2, i)) unitBx = SIN(angle) unitBy = -COS(angle) ! unitB points outward on the n1-n2 side (away from ! the n3-n4 side). delVx = v(1, n1) * fPFlt(1, 1, 1, 4, i) + v(2, n1) * fPFlt(2, 1, 1, 4, i) & & + v(1, n2) * fPFlt(1, 1, 2, 4, i) + v(2, n2) * fPFlt(2, 1, 2, 4, i) & & - v(1, n3) * fPFlt(1, 1, 2, 4, i) - v(2, n3) * fPFlt(2, 1, 2, 4, i) & & - v(1, n4) * fPFlt(1, 1, 1, 4, i) - v(2, n4) * fPFlt(2, 1, 1, 4, i) delVy = v(1, n1) * fPFlt(1, 2, 1, 4, i) + v(2, n1) * fPFlt(2, 2, 1, 4, i) & & + v(1, n2) * fPFlt(1, 2, 2, 4, i) + v(2, n2) * fPFlt(2, 2, 2, 4, i) & & - v(1, n3) * fPFlt(1, 2, 2, 4, i) - v(2, n3) * fPFlt(2, 2, 2, 4, i) & & - v(1, n4) * fPFlt(1, 2, 1, 4, i) - v(2, n4) * fPFlt(2, 2, 1, 4, i) ! delVx and delVy are the velocities of the n1-n2 side ! relative to the n3-n4 side. spread = delVx * unitBx + delVy * unitBy delTau = constr * spread tlan = 0.5D0 * (tLNode(n1) + tLNode(n2)) zman = 0.5D0 * (zMNode(n1) + zMNode(n2)) IF ((tlan <= 0.0D0).OR.(zTranF(2, i) <= 0.0D0)) THEN ! Crust alone resists convergence: dPMax = -2.0D0 * deltau / zTranF(1, i) dDPNdZ = dpmax / zTranF(1, i) ELSE ! Mantle lithosphere helps to resist convergence: dDPNdZ = -deltau / & & (0.50D0 * zTranF(1, i)**2 + zTranF(2, i) * zman + & & 0.50D0 * zTranF(2, i)**2) END IF ! dDPNdZ is the gradient of excess normal pressure (in ! excess of vertical pressure) with depth on this fault; ! check that it lies within frictional limits of blocks: q = 0.250D0 * (dQdTdA(n1) + dQdTdA(n2) + & & dQdTdA(n3) + dQdTdA(n4)) tTrans = tSurf + zTranF(1, i) * q / conduc(1) - & & zTranF(1, i)**2 * radio(1) / (2.0D0 * conduc(1)) tMeanC = (tSurf + tTrans) / 2.0D0 rhoC = rhoBar(1) * (1.0D0 - alphaT(1) * tMeanC) dLEPdC = gMean * (rhoC - rhoH2O * t_Biot) thrust = dLEPdC * cGamma normal = dLEPdC / cGamma dDPNdZ = MAX(dDPNdZ, normal - dLEPdC) dDPNdZ = MIN(dDPNdZ, thrust - dLEPdC) ELSE ! Different logic will be used; this parameter is not ! really needed. Zero it just to be careful. dDPNdZ = 0.0D0 END IF DO 90 m = 1, 7 ! elevation: elevat = elev(n1) * fPhi(1, m) + elev(n2) * fPhi(2, m) ! heat flow: q = dQdTdA(n1) * fPhi(1, m) + dQdTdA(n2) * fPhi(2, m) ! crustal thickness: crust = zMNode(n1) * fPhi(1, m) + zMNode(n2) * fPhi(2, m) ! mantle lithosphere thickness: mantle = tLNode(n1) * fPhi(1, m) + tLNode(n2) * fPhi(2, m) mantle = MAX(mantle, 0.0D0) ! Moho temperature: tMoho = tSurf + crust * q / conduc(1) - & & crust**2 * radio(1) / (2.0D0 * conduc(1)) ! Temperature at base of plate: tAsth = tMoho + mantle * (q - crust * radio(1)) / conduc(2) - & & mantle**2 * radio(2) / (2.0D0 * conduc(2)) ! mean temperatures: tMean(1) = (tSurf + tMoho) / 2.0D0 tMean(2) = (tMoho + tAsth) / 2.0D0 ! mean densities: rho(1) = rhoBar(1) * (1.0D0 - alphaT(1) * tMean(1)) rho(2) = rhoBar(2) * (1.0D0 - alphaT(2) * tMean(2)) ! derivitives of lithostatic effective pressure wrt depth dLEPdZ(1) = gMean * (rho(1) - rhoH2O * t_Biot) ePMoho = dLEPdZ(1) * crust dLEPdZ(2) = gMean * (rho(2) - rhoH2O * t_Biot) ! "angle" is the fault strike, in radians cclkws from +X. !CCCC angle = fArg(1, i) * fPhi(1, m) + fArg(2, i) * fPhi(2, m) !CCCC Line above was replaced due to cycle-shift problem! angle = Chord(fArg(1, i), fPhi(2, m), fArg(2, i)) ! unitA is a unit vector along the fault, from n1 to n2. unitAx = COS(angle) unitAy = SIN(angle) ! unitB is a perpendicular unit vector, pointing out ! toward the n4-n3 side. unitBx = -unitAy unitBy = + unitAx ! relative velocities are for n1-n2 side relative to ! the n4-n3 side: delVx = v(1, n1) * fPFlt(1, 1, 1, m, i) + v(2, n1) * fPFlt(2, 1, 1, m, i) & & + v(1, n2) * fPFlt(1, 1, 2, m, i) + v(2, n2) * fPFlt(2, 1, 2, m, i) & & - v(1, n3) * fPFlt(1, 1, 2, m, i) - v(2, n3) * fPFlt(2, 1, 2, m, i) & & - v(1, n4) * fPFlt(1, 1, 1, m, i) - v(2, n4) * fPFlt(2, 1, 1, m, i) delVy = v(1, n1) * fPFlt(1, 2, 1, m, i) + v(2, n1) * fPFlt(2, 2, 1, m, i) & & + v(1, n2) * fPFlt(1, 2, 2, m, i) + v(2, n2) * fPFlt(2, 2, 2, m, i) & & - v(1, n3) * fPFlt(1, 2, 2, m, i) - v(2, n3) * fPFlt(2, 2, 2, m, i) & & - v(1, n4) * fPFlt(1, 2, 1, m, i) - v(2, n4) * fPFlt(2, 2, 1, m, i) ! sinistral strike-slip rate component: sinist = delVx * unitAx + delVy * unitAy ! convergence rate component (in horizontal plane): close = delVx * unitBx + delVy * unitBy ! dip of the fault (from horizontal on the n1-n2 side): dip = fDip(1, i) * fPhi(1, m) + fDip(2, i) * fPhi(2, m) sloped = ABS(dip - 1.57079632679490D0) > wedge IF (.NOT.sloped) THEN ! case of a near-vertical fault: dSFdZ(1) = (dLEPdZ(1) + dDPNdZ) * fric sFMoho = dSFdZ(1) * crust dSFdZ(2) = (dLEPdZ(2) + dDPNdZ) * fric slip = ABS(sinist) locked = .FALSE. ELSE ! case of a shallow-dipping fault: ! vUpDip is the up-dip velocity component, in the ! fault plane, of the block on the n1-n3 side: vUpDip = close / COS(dip) ! "rake" angle is measured counterclockwise in ! fault plane from horizontal & parallel to "angle": rake = ATan2F(vUpDip, sinist) ! derivitive of effective normal pressure ! with respect to shear traction on fault: dEPdST = TAN(dip) * SIN(rake) ! (Notice that when sense of dip reverses, sign ! change caused by TAN(dip) is cancelled by sign ! change caused by SIN(rake).) ! According to theory, the equation to solve is: ! d(shear traction)/dz = ! fric*(dLEPdZ+dEPdST*d(shear_traction)/dz) ! This may have a physical solution (one with ! positive shear traction). If not, the fault ! is locked. locked = ((fric * dEPdST) >= 1.00D0) IF (locked) THEN dSFdZ(1) = huge dSFdZ(2) = huge ELSE dSFdZ(1) = fric * dLEPdZ(1) / (1.00D0 - fric * dEPdST) sFMoho = dSFdZ(1) * crust dSFdZ(2) = fric * dLEPdZ(2) / (1.00D0 - fric * dEPdST) END IF slip = SQRT(sinist**2 + vUpDip**2) END IF slip = MAX(slip, 1.0D8 * tiny) ! Locate plastic/creep transition(s) ! by iterated halving of domain: IF (mantle > 0.0D0) THEN limit = 2 ELSE limit = 1 zTrans(2) = 0.0D0 sheart(2) = 0.0D0 END IF DO 60 layer = 1, limit topZ = 0.0D0 IF (layer == 1) THEN baseZ = crust sf0 = 0.0D0 t0 = tSurf q0 = q z0 = 0.0D0 ELSE baseZ = mantle sf0 = sFMoho t0 = tMoho q0 = q - crust * radio(1) z0 = crust END IF DO 50 kIter = 1, 15 z = 0.50D0 * (topZ + baseZ) zAbs = z + z0 shearf = z * dSFdZ(layer) + sf0 shearp = MIN(shearf, t_dCreep(layer)) t = t0 + q0 * z / conduc(layer) - (radio(layer) / & & (2.0D0 * conduc(layer))) * z**2 IF (zAbs <= (15.0D0 * oneKm)) THEN t90pc = 0.50D0 * zAbs ELSE IF (zAbs < (45.0D0 * oneKm)) THEN t90pc = (405.0D0 / 8.0D0) * oneKm + & & (-7.0D0) * zAbs + & & (13.0D0 / 40.0D0) * oneKm * (zAbs / oneKm)**2 + & & (-1.0D0 / 300.0D0) * oneKm * (zAbs / oneKm)**3 ELSE t90pc = 2.0D0 * zAbs END IF ! See Turcotte et al (1980) JGR, 85, B11, 6224-6230. strain = slip / t90pc shearc = t_aCreep(layer) * (strain**t_eCreep) * & & EXP((t_bCreep(layer) + t_cCreep(layer) * z) / t) IF (shearc < shearp) THEN baseZ = z ELSE topZ = z END IF 50 CONTINUE zTrans(layer) = 0.50D0 * (topZ + baseZ) sheart(layer) = zTrans(layer) * dSFdZ(layer) + sf0 60 CONTINUE ! plastic part of vertical integral(s) of traction: ! (A) crust: IF (sheart(1) <= t_dCreep(1)) THEN vitdz = 0.50D0 * sheart(1) * zTrans(1) ELSE zp = zTrans(1) * t_dCreep(1) / sheart(1) vitdz = t_dCreep(1) * (zTrans(1) - 0.50D0 * zp) END IF ! (B) mantle lithosphere: IF ((mantle > 0.).AND.(sheart(2) > sfmoho)) THEN IF (sheart(2) <= t_dCreep(2)) THEN vitdz = vitdz + 0.50D0 * (sfmoho + sheart(2)) * zTrans(2) ELSE zp = zTrans(2) * (t_dCreep(2) - sfmoho) / & & (sheart(2) - sfmoho) zp = MAX(zp, 0.) vitdz = vitdz + 0.50D0 * (sfmoho + sheart(2)) * zp + & & t_dCreep(2) * (zTrans(2) - zp) END IF END IF ! Add creep part(s) of integral, using parabolic rule: sum = 0.0D0 DO 80 layer = 1, limit IF (layer == 1) THEN thick = crust t0 = tSurf q0 = q zAbs = 0.0D0 ELSE thick = mantle t0 = tMoho q0 = q - crust * radio(1) zAbs = crust END IF dz = (thick - zTrans(layer)) / nStep oldsc = sheart(layer) oldsc = MIN(oldsc, t_dCreep(layer)) z0 = zTrans(layer) DO 70 j = 1, nStep zhalf = z0 + 0.50D0 * dz zfull = z0 + dz azhalf = zhalf + zAbs azfull = zfull + zAbs thalf = t0 + q0 * zhalf / conduc(layer) - & & (radio(layer) / & & (2.0D0 * conduc(layer))) * zhalf**2 tfull = t0 + q0 * zfull / conduc(layer) - & & (radio(layer) / & & (2.0D0 * conduc(layer))) * zfull**2 IF (azhalf <= (15.0D0 * oneKm)) THEN whalf = 0.50D0 * azhalf ELSE IF (azhalf < (45.0D0 * oneKm)) THEN whalf = (405.0D0 / 8.0D0) * oneKm + & & (-7.0D0) * azhalf + & & (13.0D0 / 40.0D0) * oneKm * (azhalf / oneKm)**2 + & & (-1.0D0 / 300.0D0) * oneKm * (azhalf / oneKm)**3 ELSE whalf = 2.0D0 * azhalf END IF IF (azfull <= (15.0D0 * oneKm)) THEN wfull = 0.50D0 * azfull ELSE IF (azfull < (45.0D0 * oneKm)) THEN wfull = (405.0D0 / 8.0D0) * oneKm + & & (-7.0D0) * azfull + & & (13.0D0 / 40.0D0) * oneKm * (azfull / oneKm)**2 + & & (-1.0D0 / 300.0D0) * oneKm * (azfull / oneKm)**3 ELSE wfull = 2.0D0 * azhalf END IF ! See Turcotte et al (1980) JGR, 85, B11, 6224-6230. ehalf = slip / whalf efull = slip / wfull schalf = t_aCreep(layer) * (ehalf**t_eCreep) * & & EXP((t_bCreep(layer) + t_cCreep(layer) * zhalf) & & / thalf) schalf = MIN(schalf, t_dCreep(layer)) scfull = t_aCreep(layer) * (efull**t_eCreep) * & & EXP((t_bCreep(layer) + t_cCreep(layer) * zfull) & & / tfull) scfull = MIN(scfull, t_dCreep(layer)) sum = sum + dz * (0.1666667D0 * oldsc + & & 0.6666667D0 * schalf + & & 0.1666666D0 * scfull) z0 = zfull oldsc = scfull 70 CONTINUE 80 CONTINUE vitdz = vitdz + sum ! Limit shear traction on subduction zones only: dippy = MIN(dip, 3.14159265358979D0 - dip) IF (dippy <= slide) THEN IF (elevat < 0.0D0) THEN ! apply oceanic subduction zone limit: vitdz = MIN(vitdz, tauMax(1) * SIN(dip)) ELSE ! apply continental subduction zone limit: vitdz = MIN(vitdz, tauMax(2) * SIN(dip)) END IF END IF dpt1 = (1.D0 * vitdz) / slip vIMuDZ = MIN(dpt1, 1.D38) fIMuDZ(m, i) = MIN(vimudz, fMuMax * (crust + mantle)) ! Dipping, oblique-slip integration ! points are also characterized ! by fC and fTStar: IF (sloped) THEN ts = sinist * fIMuDZ(m, i) tu = vUpDip * fIMuDZ(m, i) IF (locked) THEN fC(1, 1, m, i) = fIMuDZ(m, i) fC(1, 2, m, i) = 0.0D0 fC(2, 1, m, i) = 0.0D0 fC(2, 2, m, i) = fIMuDZ(m, i) ELSE sinr = SIN(rake) cosr = COS(rake) tand = TAN(dip) ! *** IMPORTANT NOTE: *** ! The following 7 statements are *not* the ! result of theory, but a tactical choice ! which attempts to compromise between ! stability of the linear system, stability ! of the iteration, and efficiency. ! They may be changed if the program does ! no converge satisfactorily! tune = 2.0D0 fC(1, 1, m, i) = fIMuDZ(m, i) * & & (1.0D0 - tune * sinr * cosr**2 * tand) fC(1, 2, m, i) = fIMuDZ(m, i) * & & (tune * cosr**3 * tand) fC(2, 1, m, i) = fIMuDZ(m, i) * & & (-tune * sinr**2 * cosr * tand) fC(2, 2, m, i) = fIMuDZ(m, i) * & & (1.0D0 + tune * sinr * cosr**2 * tand) ! (Often, fC(1,2) is the biggest term. ! In some cases, diagonals become negative. ! For stability, be sure that the fC ! matrix remains positive-definite: fC(1, 1, m, i) = MAX(fC(1, 1, m, i), ABS(fC(1, 2, m, i))) fC(2, 2, m, i) = MAX(fC(2, 2, m, i), ABS(fC(1, 2, m, i))) END IF fTStar(1, m, i) = ts - fC(1, 1, m, i) * sinist - & & fC(1, 2, m, i) * vUpDip fTStar(2, m, i) = tu - fC(2, 1, m, i) * sinist - & & fC(2, 2, m, i) * vUpDip END IF ! Provide interesting diagnostic data at midpoints only: IF (m == 4) THEN fSlips(i) = (.NOT.locked).AND. & & (fIMuDZ(m, i) < (0.99D0 * fMuMax * (crust + mantle))) zTranF(1, i) = zTrans(1) fPeakS(1, i) = MIN(sheart(1), t_dCreep(1)) zTranF(2, i) = zTrans(2) fPeakS(2, i) = MIN(sheart(2), t_dCreep(2)) END IF 90 CONTINUE 100 CONTINUE RETURN END SUBROUTINE Mohr SUBROUTINE Next (i, j, mxEl, mxFEl, nFl, nodeF, nodes, numEl, & ! input & kFault, kEle) ! output ! Determine whether there are more elements adjacent to side #j of ! triangular continuum element #i. ! j = 1 means the side opposite node # nodes(1, i). ! j = 2 means the side opposite node # nodes(2, i). ! j = 3 means the side opposite node # nodes(3, i). ! If a fault element is adjacent, its number is kFault; otherwise, ! kFault is set to zero. ! If another triangular continuum element is adjacent (even across ! fault element kFault !) then its number is kEle; otherwise, kEle = 0. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - INTEGER, INTENT(IN) :: i, j, mxEl, mxFEl, nFl, nodeF, nodes, numEl ! input INTEGER, INTENT(OUT) :: kFault, kEle ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - INTEGER k, khigh, klow, l, m1, m2, m3, m4, n1, n2 LOGICAL foundF DIMENSION nodeF(4, mxFEl), nodes(3, mxEl) ! Two node numbers along the side of interest, counterclockwise: n1 = nodes(MOD(j, 3) + 1, i) n2 = nodes(MOD(j + 1, 3) + 1, i) ! Check for adjacent fault element first: foundF = .FALSE. kFault = 0 IF (nFl > 0) THEN DO 10 k = 1, nFl m1 = nodeF(1, k) m2 = nodeF(2, k) m3 = nodeF(3, k) m4 = nodeF(4, k) IF (((m1 == n2).AND.(m2 == n1)).OR. & & ((m3 == n2).AND.(m4 == n1))) THEN foundF = .TRUE. kFault = k GO TO 11 END IF 10 CONTINUE 11 IF (.NOT.foundF) kFault = 0 ! If there was a fault, replace 2 node numbers that we search for: IF (foundF) THEN IF (m2 == n1) THEN n1 = m3 n2 = m4 ELSE n1 = m1 n2 = m2 END IF END IF END IF ! Search for adjacent triangular continuum element: kEle = 0 klow = i khigh = i ! --- Begin irregular loop, to search out nearest elements first --- 100 klow = klow - 1 IF (klow >= 1) THEN DO 110 l = 1, 3 m1 = nodes(MOD(l, 3) + 1, klow) m2 = nodes(MOD(l + 1, 3) + 1, klow) IF ((m2 == n1).AND.(m1 == n2)) THEN kEle = klow RETURN END IF 110 CONTINUE END IF khigh = khigh + 1 IF (khigh <= numEl) THEN DO 120 l = 1, 3 m1 = nodes(MOD(l, 3) + 1, khigh) m2 = nodes(MOD(l + 1, 3) + 1, khigh) IF ((m2 == n1).AND.(m1 == n2)) THEN kEle = khigh RETURN END IF 120 CONTINUE END IF IF ((klow > 1).OR.(khigh < numEl)) GO TO 100 RETURN END SUBROUTINE Next SUBROUTINE OldVel (iUnitT, iUnitV, mxNode, numNod, & ! INTENT(IN) & haveNV, title1, title2, title3, v) ! INTENT(OUT) ! READ old velocity solution from unit iUnitV, or else set LOGICAL ! variable 'haveNV" to .FALSE. ! WRITEs 3 title lines to iUnitT. IMPLICIT NONE INTEGER, INTENT(IN) :: iUnitT, iUnitV, mxNode, numNod LOGICAL, INTENT(OUT) :: haveNV CHARACTER*80, INTENT(OUT) :: title1, title2, title3 REAL*8, DIMENSION(2, mxNode), INTENT(OUT) :: v INTEGER :: i, j READ (iUnitV, '(A80)', END = 100, ERR = 100) title1 READ (iUnitV, '(A80)', END = 100, ERR = 100) title2 READ (iUnitV, '(A80)', END = 100, ERR = 100) title3 READ (iUnitV, * , END = 100, ERR = 100) ((v(j, i), j = 1, 2), i = 1, numNod) haveNV = .TRUE. WRITE (iUnitT, 50) iUnitV, title1, title2, title3 50 FORMAT (/ /' Velocity solution was', & & ' READ from unit', I3, '; titles were:' / 3(/ ' ', A80) ) GO TO 900 ! ------------------(This section executed only if READ fails)--------- 100 WRITE (iUnitT, 110) iUnitV 110 FORMAT (/ /' NO FURTHER VELOCITY SOLUTIONS FOUND ON UNIT', I3) haveNV = .FALSE. ! --------------------------------------------------------------------- 900 WRITE (iUnitT, 999) 999 FORMAT (' -------------------------------------------------------------------------------') END SUBROUTINE OldVel SUBROUTINE OnArc (s, theta1, phi1, theta2, phi2, & ! INTENT(IN) & theta, phi) ! INTENT(OUT) ! Computes coordinates (theta, phi) = (colatitude, longitude) in radians ! for a point which lies on the great circle arc from (theta1, phi1) to (theta2, phi2). ! Input parameter s expresses the fractional distance, ! so input s = 0.0D0 returns (theta1, phi1) and ! input s = 1.0D0 returns (theta2, phi2). IMPLICIT NONE REAL*8, INTENT(IN) :: s, theta1, phi1, theta2, phi2 REAL*8, INTENT(OUT) :: theta, phi REAL*8 :: comple, equat, size REAL*8 :: tvec(3), uvec(3), uvec1(3), uvec2(3) ! These are all unit vectors (in the unit sphere) ! in a Cartesian coordinate system with ! X = (0E, 0N), Y = (90E, 0N), Z = North pole. uvec1(1) = SIN(theta1) * COS(phi1) uvec1(2) = SIN(theta1) * SIN(phi1) uvec1(3) = COS(theta1) uvec2(1) = SIN(theta2) * COS(phi2) uvec2(2) = SIN(theta2) * SIN(phi2) uvec2(3) = COS(theta2) comple = 1.0D0 - s tvec(1) = comple * uvec1(1) + s * uvec2(1) tvec(2) = comple * uvec1(2) + s * uvec2(2) tvec(3) = comple * uvec1(3) + s * uvec2(3) size = SQRT(tvec(1)**2 + tvec(2)**2 + tvec(3)**2) uvec(1) = tvec(1) / size uvec(2) = tvec(2) / size uvec(3) = tvec(3) / size equat = SQRT(uvec(1)**2 + uvec(2)**2) theta = ATAN2(equat, uvec(3)) phi = ATan2F(uvec(2), uvec(1)) END SUBROUTINE OnArc SUBROUTINE Pause() ! This routine is called after any error message, ! and also after successful program completion, ! in order to let the user have time to read the message ! when the job is run interactively. ! To run the job as a batch job (e.g., from a Windows .bat file), ! comment out the READ so that this routine does not stop ! the flow of the program. IMPLICIT NONE WRITE(*, "(' Press [Enter]...')") ! ------------------------------------- !CCCC Following line could be commented out so that jobs won't pause !CCCC on batch-processing computer systems. !CCCC Remove the !CCCC for use on interactive computer systems! !CCCC READ(*, *) ! ------------------------------------- END SUBROUTINE Pause SUBROUTINE Prince (e11, e22, e12, & ! input & e1, e2, u1x, u1y, u2x, u2y) ! output ! Find principal values (e1, e2) of the symmetric 2x2 tensor e11 e12 ! e12 e22 ! and also the associated eigenvectors #1=(u1x, u1y); #2=(u2x, u2y). ! The convention is that e1 <= e2. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - REAL*8, INTENT(IN) :: e11, e22, e12 ! input REAL*8, INTENT(OUT) :: e1, e2, u1x, u1y, u2x, u2y ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - REAL*8 c, r, scale, test, theta r = SQRT(((e11 - e22) * 0.5D0)**2 + e12**2) c = (e11 + e22) * 0.5D0 e1 = c - r e2 = c + r scale = MAX(ABS(e1), ABS(e2)) test = 0.01D0 * scale IF ((ABS(e12) > test).OR.(ABS(e11 - e1) > test)) THEN theta = ATan2F(e11 - e1, -e12) ELSE theta = ATan2F(e12, e1 - e22) END IF u1x = COS(theta) u1y = SIN(theta) u2x = u1y u2y = -u1x RETURN END SUBROUTINE Prince SUBROUTINE Projec (iEle, mxEl, mxNode, nodes, & & s1, s2, s3, & & xNode, yNode, v, & ! INTENT(IN) & vTheta, vPhi) ! INTENT(OUT) ! Computes horizontal velocity (vTheta,vPhi) at point (s1, s2, s3) ! in triangular element iEle. IMPLICIT NONE INTEGER, INTENT(IN) :: iEle, mxEl, mxNode INTEGER, DIMENSION(3, mxEl), INTENT(IN) :: nodes REAL*8, INTENT(IN) :: s1, s2, s3 REAL*8, DIMENSION(mxNode), INTENT(IN) :: xNode, yNode REAL*8, DIMENSION(2, mxNode), INTENT(IN) :: v REAL*8, INTENT(OUT) :: vTheta, vPhi INTEGER :: k REAL*8 :: growth, size REAL*8 :: rn(3, 3), rP(3), rS(3), & & uTheta(3, 3), uPhi(3, 3), uT(3), uP(3), & & vCart(3, 3), vP(3), vS(3) DO 10 k = 1, 3 ! Find Cartesian radius vector for each node: rn(1, k) = SIN(xNode(nodes(k, iEle))) * COS(yNode(nodes(k, iEle))) rn(2, k) = SIN(xNode(nodes(k, iEle))) * SIN(yNode(nodes(k, iEle))) rn(3, k) = COS(xNode(nodes(k, iEle))) ! Find Theta-pointing Cartesian unit vector at each node: uTheta(1, k) = COS(xNode(nodes(k, iEle))) * COS(yNode(nodes(k, iEle))) uTheta(2, k) = COS(xNode(nodes(k, iEle))) * SIN(yNode(nodes(k, iEle))) uTheta(3, k) = -SIN(xNode(nodes(k, iEle))) ! Find Phi-pointing Cartesian unit vector at each node: uPhi(1, k) = -SIN(yNode(nodes(k, iEle))) uPhi(2, k) = COS(yNode(nodes(k, iEle))) uPhi(3, k) = 0.0D0 ! Create Cartesian velocity at each node: vCart(1, k) = v(1, nodes(k, iEle)) * uTheta(1, k) + & & v(2, nodes(k, iEle)) * uPhi(1, k) vCart(2, k) = v(1, nodes(k, iEle)) * uTheta(2, k) + & & v(2, nodes(k, iEle)) * uPhi(2, k) vCart(3, k) = v(1, nodes(k, iEle)) * uTheta(3, k) + & & v(2, nodes(k, iEle)) * uPhi(3, k) 10 CONTINUE ! Interpolate velocity in plane triangle: vP(1) = s1 * vCart(1, 1) + s2 * vCart(1, 2) + s3 * vCart(1, 3) vP(2) = s1 * vCart(2, 1) + s2 * vCart(2, 2) + s3 * vCart(2, 3) vP(3) = s1 * vCart(3, 1) + s2 * vCart(3, 2) + s3 * vCart(3, 3) ! Interpolate position in plane triangle: rP(1) = s1 * rn(1, 1) + s2 * rn(1, 2) + s3 * rn(1, 3) rP(2) = s1 * rn(2, 1) + s2 * rn(2, 2) + s3 * rn(2, 3) rP(3) = s1 * rn(3, 1) + s2 * rn(3, 2) + s3 * rn(3, 3) size = SQRT(rP(1)**2 + rP(2)**2 + rP(3)**2) growth = 1.0D0 / size ! Project position to sphere: rS(1) = rP(1) * growth rS(2) = rP(2) * growth rS(3) = rP(3) * growth ! Project velocity to velocity on sphere: vS(1) = vP(1) * growth vS(2) = vP(2) * growth vS(3) = vP(3) * growth ! Find Theta and Phi unit vectors at this point: uP(1) = -rS(2) uP(2) = rS(1) uP(3) = 0.0D0 CALL Unit (uP) CALL DCross (uP, rS, uT) ! Find horizontal components of velocity: vTheta = vS(1) * uT(1) + vS(2) * uT(2) + vS(3) * uT(3) vPhi = vS(1) * uP(1) + vS(2) * uP(2) + vS(3) * uP(3) END SUBROUTINE Projec SUBROUTINE PutNet (iUnitO, & ! INTENT(IN) & brief, dQdTdA, elev, fDip, & & mxEl, mxFEl, mxNode, n1000, & & nFakeN, nFl, nodeF, nodes, & & nRealN, numEl, numNod, offset, & & title1, tLNode, xNode, yNode, zMNode) ! WRITEs finite element grid to unit "iUnitO". USE DSphere ! in Peter Bird's file DSphere.f90 IMPLICIT NONE INTEGER, INTENT(IN) :: iUnitO LOGICAL, INTENT(IN) :: brief REAL*8, DIMENSION(mxNode), INTENT(IN) :: dQdTdA, elev REAL*8, DIMENSION(2, mxFEl), INTENT(IN) :: fDip INTEGER, INTENT(IN) :: mxEl, mxFEl, mxNode, n1000, nFakeN, nFl INTEGER, DIMENSION(4, mxFEl), INTENT(IN) :: nodeF INTEGER, DIMENSION(3, mxEl), INTENT(IN) :: nodes INTEGER, INTENT(IN) :: nRealN, numEl, numNod REAL*8, DIMENSION(mxFEl), INTENT(IN) :: offset CHARACTER*80, INTENT(IN) :: title1 REAL*8, DIMENSION(mxNode), INTENT(IN) :: tLNode, xNode, yNode, zMNode INTEGER :: i, iPrint, k, np(4) REAL*8 :: dips(2), pLat, pLon WRITE (iUnitO, 1) title1 1 FORMAT (A80) WRITE (iUnitO, 2) numNod, nRealN, nFakeN, n1000, brief 2 FORMAT(4I8, L8, ' (numNod, nRealN, nFakeN, N1000, BRIEF)') DO 100 i = 1, numNod IF (i <= nRealN) THEN iPrint = i ELSE iPrint = n1000 + (i - nRealN) END IF pLat = 90.0D0 - xNode(i) * degrees_per_radian pLon = yNode(i) * degrees_per_radian WRITE (iUnitO, 91) i, pLon, pLat, elev(i), dQdTdA(i), zMNode(i), tLNode(i) 91 FORMAT (I8, 2F11.5, 4ES10.2) 100 CONTINUE WRITE (iUnitO, 110) numEl 110 FORMAT (I10,' (numEl = NUMBER OF TRIANGULAR CONTINUUM ELEMENTS)') DO 200 i = 1, numEl DO 150 k = 1, 3 IF (nodes(k, i) <= nRealN) THEN np(k) = nodes(k, i) ELSE np(k) = n1000 + (nodes(k, i) - nRealN) END IF 150 CONTINUE WRITE (iUnitO, 160) i, (np(k), k = 1, 3) 160 FORMAT (I8, 3I8) 200 CONTINUE WRITE (iUnitO, 210) nFl 210 FORMAT (I10,' (nFl = NUMBER OF CURVILINEAR FAULT ELEMENTS)') DO 300 i = 1, nFl DO 220 k = 1, 4 IF (nodeF(k, i) <= nRealN) THEN np(k) = nodeF(k, i) ELSE np(k) = n1000 + (nodeF(k, i) - nRealN) END IF 220 CONTINUE DO 230 k = 1, 2 dips(k) = fDip(k, i) dips(k) = dips(k) * degrees_per_radian IF (dips(k) > 90.01D0) dips(k) = dips(k) - 180.0D0 230 CONTINUE WRITE (iUnitO, 250) i, (np(k), k = 1, 4), (dips(k), k = 1, 2), offset(i) 250 FORMAT (I8, 4I6, 2F5.0, ES10.2) 300 CONTINUE END SUBROUTINE PutNet SUBROUTINE ReadN (iUnitP, iUnitT, n, & ! input & vector) ! output ! A utility routine designed to permit -Faults- input files ! to also be used by -Plates-, which expects more numbers ! in some records. ! This routine attempts to READ "n" floating-point values ! (using * FORMAT) from the next record on device iUnitP. ! If anything goes wrong, the missing values are set to zero. ! Since June 2005, this feature also allows the reading of ! both old-format (OrbData) .feg files (with 4 nodal data), ! and new-format (OrbData5) .feg files (with 6 nodal data), ! by program -Shells-. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - INTEGER, INTENT(IN) :: iUnitP, iUnitT, n ! input REAL*8, INTENT(OUT) :: vector ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - CHARACTER*1 c CHARACTER*132 line INTEGER i, number LOGICAL anyIn, dotted, expon, signed DIMENSION vector(n) line = ' ' READ (iUnitP, 1) line 1 FORMAT (A) number = 0 anyIn = .FALSE. expon = .FALSE. signed = .FALSE. dotted = .FALSE. DO 10 i = 1, 132 c = line(i:i) IF ((c == ' ').OR.(c == ',').OR.(c == '/')) THEN signed = .FALSE. expon = .FALSE. dotted = .FALSE. IF (anyIn) THEN number = number + 1 anyIn = .FALSE. END IF ELSE IF ((c == '+').OR.(c == '-')) THEN IF (signed) THEN GO TO 50 ELSE signed = .TRUE. END IF ELSE IF ((c == 'E').OR.(c == 'D').OR. & & (c == 'e').OR.(c == 'd')) THEN IF (expon) THEN GO TO 50 ELSE expon = .TRUE. signed = .FALSE. dotted = .TRUE. END IF ELSE IF (c == '.') THEN IF (dotted) THEN GO TO 50 ELSE dotted = .TRUE. END IF ELSE IF ((c == '0').OR.(c == '1').OR.(c == '2').OR. & & (c == '3').OR.(c == '4').OR.(c == '5').OR. & & (c == '6').OR.(c == '7').OR.(c == '8').OR. & & (c == '9')) THEN signed = .TRUE. anyIn = .TRUE. ELSE GO TO 50 END IF 10 CONTINUE IF (anyIn) number = number + 1 50 IF (number == 0) THEN WRITE (iUnitT, 91) n, line 91 FORMAT (/' ERR0R: A LINE OF PARAMETER INPUT WHICH', & & ' WAS SUPPOSED TO CONTAIN 1-', I2, ' NUMBERS'/ & & ' COULD NOT BE INTERPRETED. LINE FOLLOWS:'/ & & ' ', A80) CALL Pause() STOP ELSE number = MIN(number, n) BACKSPACE iUnitP READ (iUnitP, * ) (vector(i), i = 1, number) IF (number < n) THEN DO 99 i = number + 1, n vector(i) = 0.0D0 99 CONTINUE END IF END IF RETURN END SUBROUTINE ReadN SUBROUTINE ReadPm (iUnit7, iUnit8, names , numPlt, offMax, & ! input & alphaT, conduc, & ! output & d_fFric, d_cFric, d_Biot, d_Byerly, d_aCreep, d_bCreep, d_cCreep, d_dCreep, d_eCreep, & & everyP, gMean , & & gradie, iConve, iPVRef, & & maxItr, okDelV, okToQt, oneKm, radio, & & radius, refStr, rhoAst, rhoBar, rhoH2O, & & tAdiab, tauMax, temLim, title3, & & trHMax, tSurf, vTimes, zBAsth) ! Reads strategic and tactical input parameters from device iUnit7, ! and echoes them on device iUnit8 with annotations. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - INTEGER, INTENT(IN) :: iUnit7, iUnit8 ! input CHARACTER*2, INTENT(IN) :: names ! input INTEGER, INTENT(IN) :: numPlt ! input REAL*8, INTENT(IN) :: offMax ! input REAL*8, INTENT(OUT) :: alphaT, conduc ! output REAL*8, INTENT(OUT) :: d_fFric, d_cFric, d_Biot, d_Byerly, d_aCreep, d_bCreep, d_cCreep, d_dCreep, d_eCreep ! output LOGICAL, INTENT(OUT) :: everyP ! output REAL*8, INTENT(OUT) :: gMean , gradie ! output INTEGER, INTENT(OUT) :: iConve, iPVRef, maxItr ! output REAL*8, INTENT(OUT) :: okDelV, okToQt, oneKm, radio, radius, refStr, & ! output & rhoAst, rhoBar, rhoH2O, tAdiab, tauMax, temLim ! output CHARACTER*80, INTENT(OUT) :: title3 ! output REAL*8, INTENT(OUT) :: trHMax, tSurf, vTimes, zBAsth ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - CHARACTER*2 pltRef INTEGER i, ios REAL*8 tempV, vector DIMENSION alphaT(2), conduc(2), & & d_aCreep(2), d_bCreep(2), d_cCreep(2), d_dCreep(2), & & names(numplt), radio(2), & & rhoBar(2), tauMax(2), temLim(2), tempv(2), vector(2) WRITE (*, 1) iUnit7 1 FORMAT(//' Attempting to read input parameter file', & & ' from unit ',I3/) title3 = ' ' READ (iUnit7, 2, IOSTAT = ios) title3 IF (ios /= 0) THEN WRITE(*, "(' ERR','OR: File not found, or file is empty,' & & /' or file is too short.')") CALL Pause() STOP END IF 2 FORMAT (A80) WRITE (iUnit8, 3) title3 3 FORMAT (/' [OMIT from iXXX.in] Title of parameter set ='/' ',A80) WRITE (iUnit8, 4) 4 FORMAT (' [OMIT from iXXX.in]' & & /' [OMIT from iXXX.in]', & & ' **************************************************' & & /' [OMIT from iXXX.in]', & & ' It is the user''s responsibility to input all of the' & & /' [OMIT from iXXX.in]', & & ' following numerical quantities in consistent units,' & & /' [OMIT from iXXX.in]', & & ' such as Systeme International (SI) or cm-g-s (cgs).' & & /' [OMIT from iXXX.in]', & & ' Note that time unit must be the second (hard-coded).' & & /' [OMIT from iXXX.in]', & & ' **************************************************' & & /' [OMIT from iXXX.in]' & & /' [OMIT from iXXX.in]', & & ' ========== Strategic parameters (define the real', & & '-Earth problem) ======' & & /' [OMIT from iXXX.in]') READ (iUnit7, * , IOSTAT = ios) d_fFric IF (ios /= 0) THEN WRITE(*, "(' ERR','OR: File not found, or file is empty,' & & /' or file is too short.')") CALL Pause() STOP END IF WRITE (iUnit8, 20) d_fFric 20 FORMAT (' ',11X,F10.3,' fFric = coefficient of friction on faults') IF (d_fFric < 0.0D0) THEN WRITE (*, 21) WRITE (iUnit8, 21) 21 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' negative fault friction fFric is not physical.') CALL Pause() STOP END IF READ (iUnit7, * ) d_cFric WRITE (iUnit8, 30) d_cFric 30 FORMAT (' ',11X,F10.3,' cFric = coefficient of friction within blocks') IF (d_cFric <= 0.0D0) THEN WRITE (*, 31) WRITE (iUnit8, 31) 31 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' continuum friction cFric must be positive.') CALL Pause() STOP END IF READ (iUnit7, * ) d_Biot WRITE (iUnit8, 40) d_Biot 40 FORMAT (' ',11X,F10.4,' Biot = effective-pressure coefficient', & & ' of Biot: 0. (dry) to 1. (wet)') IF ((d_Biot < 0.0D0).OR.(d_Biot > 1.0D0)) THEN WRITE (*, 41) WRITE (iUnit8, 41) 41 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' Biot coefficient must be in range 0.0 to 1.0.') CALL Pause() STOP END IF READ (iUnit7, * ) d_Byerly IF (offMax > 0.0D0) THEN WRITE (iUnit8, 43) d_Byerly 43 FORMAT (' ',11X,F10.4,' Byerlee coefficient (0. to 0.999) ='/ & & 11X,' fractional friction reduction on master fault'/ & & 11X,' (Other faults have less reduction, in proportion to'/ & & 11X,' their total past offsets)') IF ((d_Byerly < 0.0D0).OR.(d_Byerly > 1.0D0)) THEN WRITE (*, 44) WRITE (iUnit8, 44) 44 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' Byerlee coefficient must be in range 0.0 to 1.0.') CALL Pause() STOP END IF ELSE WRITE (iUnit8, 46) d_Byerly 46 FORMAT (' ',11X,F10.4,' Byerlee coefficient (not used in', & & ' this run, as all fault offsets are zero).') END IF CALL ReadN (iUnit7, iUnit8, 2, & ! input & d_aCreep) ! output WRITE (iUnit8, 50) d_aCreep(1), d_aCreep(2) 50 FORMAT (' ',1P, E10.2,' ',E10.2,' aCreep = A for creep = ', & & 'pre-exponential shear', & & ' stress constant for creep. (crust/mantle)') IF ((d_aCreep(1) <= 0.0D0).OR.(d_aCreep(2) <= 0.0D0)) THEN WRITE (*, 51) WRITE (iUnit8, 51) 51 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' aCreep must be positive in each layer.') CALL Pause() STOP END IF CALL ReadN (iUnit7, iUnit8, 2, & ! input & d_bCreep) ! output WRITE (iUnit8, 60) d_bCreep(1), d_bCreep(2) 60 FORMAT (' ',F10.0,' ',F10.0,' bCreep = B for creep =', & & ' activation_energy/R/n', & & ' (in K). (crust/mantle)') IF ((d_bCreep(1) < 0.0D0).OR.(d_bCreep(2) < 0.0D0)) THEN WRITE (*, 61) WRITE (iUnit8, 61) 61 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' Negative bCreep in either layer is unphysical.') CALL Pause() STOP END IF CALL ReadN (iUnit7, iUnit8, 2, & ! input & d_cCreep) ! output WRITE (iUnit8, 70) d_cCreep(1), d_cCreep(2) 70 FORMAT (' ',1P, E10.2,' ',E10.2,' cCreep = C for creep = rho*', & & 'g*V_star*eCreep/R (in K/m). (crust/mantle)') IF ((d_cCreep(1) < 0.0D0).OR.(d_cCreep(2) < 0.0D0)) THEN WRITE (*, 71) WRITE (iUnit8, 71) 71 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' Negative cCreep in either layer is unphysical.') CALL Pause() STOP END IF CALL ReadN (iUnit7, iUnit8, 2, & ! input & d_dCreep) ! output WRITE (iUnit8, 80) d_dCreep(1), d_dCreep(2) 80 FORMAT (' ',1P,E10.2,' ',E10.2,' dCreep = D for creep = max', & & 'imum shear stress under any conditions. (crust/mantle)') IF ((d_dCreep(1) <= 0.0D0).OR.(d_dCreep(2) <= 0.0D0)) THEN WRITE (*, 81) WRITE (iUnit8, 81) 81 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' dCreep must be positive in each layer.') CALL Pause() STOP END IF READ (iUnit7, * ) d_eCreep WRITE (iUnit8, 90) d_eCreep 90 FORMAT (' ',11X,F10.6,' eCreep = E for creep = strain-rate expo', & & 'nent for creep (1/n). (Same for crust and mantle!)') IF (d_eCreep <= 0.0D0) THEN WRITE (*, 91) WRITE (iUnit8, 91) 91 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' eCreep must be positive.') CALL Pause() STOP END IF READ (iUnit7, * ) tAdiab, gradie WRITE (iUnit8, 92) tAdiab, gradie 92 FORMAT (' ',F10.0,' ',1P,E10.2,' tAdiab, GRADIE = intercept and ' & & ,'slope of upper mantle adiabat below plate (K, K/m)') IF ((tAdiab < 0.0D0).OR.(gradie < 0.0D0)) THEN WRITE (*, 93) WRITE (iUnit8, 93) 93 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' Negative Kelvin temperature and/or' & & /' negative adiabatic gradient is/are unphysical.') CALL Pause() STOP END IF READ (iUnit7, * ) zBAsth WRITE (iUnit8, 94) zBAsth 94 FORMAT (' ',11X,1P,E10.2,' zBAsth = depth of base of', & & ' asthenosphere') IF (zBAsth <= 0.0D0) THEN WRITE (*, 95) WRITE (iUnit8, 95) 95 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' zBAsth must be positive.') CALL Pause() STOP END IF READ (iUnit7, 952) pltRef 952 FORMAT(A2) WRITE (iUnit8, 954) pltRef 954 FORMAT(' ',A2,'<==================', & & ' pltRef = plate defining velocity ', & & 'reference frame (AF, NA, EU, ...)') iPVRef = 0 DO 956 i = 1, numPlt IF (names(i) == pltRef) iPVRef = i 956 CONTINUE IF (iPVRef == 0) THEN WRITE (*, 958) (names(i), i = 1, numPlt) WRITE (iUnit8, 958) (names(i), i = 1, numPlt) 958 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' In line 13 (after zBAsth, before iConve),' & & /' in the first two columns of the line,' & & /' define the velocity reference frame by' & & /' entering one of the following plate names:' & & /' ',26(A2,1X)) CALL Pause() STOP END IF READ (iUnit7, * ) iConve WRITE (iUnit8, 96) iConve 96 FORMAT (' ',11X,I10,' iConve = code for mantle flow below the ', & & 'asthenosphere:' & &/' ','[OMIT from iXXX.in] ',' 0 = static (with respect to AF)' & &/' ','[OMIT from iXXX.in] ',' 1 = Hager and O''Connell (1979)', & & ' Model II' & &/' ','[OMIT from iXXX.in] ',' 2 = Baumgardner (1988) Figure', & & ' 7A-F, *10.' & &/' ','[OMIT from iXXX.in] ',' 3 = PB2002 (Bird, 2003)' & &/' ','[OMIT from iXXX.in] ',' 4 = PB2002 drags continents;', & & ' no ocean drag' & &/' ','[OMIT from iXXX.in] ',' 5 = drag on base of subduction', & & ' forearc only' & &/' ','[OMIT from iXXX.in] ',' 6 = sense & traction from trac', & & 'tion pole vectors' & & ) IF ((iConve < 0).OR.(iConve > 6)) THEN WRITE (*, 97) WRITE (iUnit8, 97) 97 FORMAT (/' *** ERR','OR in input parameter file: ***' & & /' iConve is not one of the values listed.') CALL Pause() STOP END IF IF (iConve > 0) THEN BACKSPACE iUnit7 CALL ReadN (iUnit7, iUnit8, 2, & ! input & tempv) ! output IF (tempv(2) >= 0) THEN vTimes = tempv(2) WRITE (iUnit8, 98) vTimes 98 FORMAT (' ',11X,F10.4,' vTimes = speed factor for con', & & 'vection model identified above') ELSE WRITE (*, 99) WRITE (iUnit8, 99) 99 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' Uninterpretable value for vTimes.') CALL Pause() STOP END IF ELSE vTimes = 1.0D0 END IF READ (iUnit7, * ) trHMax WRITE (iUnit8, 101) trHMax 101 FORMAT (' ',11X,1P,E10.2,' trHMax = limit on horizontal', & & ' tractions applied to base of plate') IF (trHMax < 0.0D0) THEN WRITE (*, 102) WRITE (iUnit8, 102) 102 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' trHMax may not be negative.') CALL Pause() STOP END IF CALL ReadN (iUnit7, iUnit8, 2, & ! input & vector) ! output tauMax(1) = vector(1) tauMax(2) = vector(2) ! Provide for old parameter files with only one tauMax: IF (tauMax(2) <= 0.0D0) tauMax(2) = tauMax(1) WRITE (iUnit8, 106) tauMax(1), tauMax(2) 106 FORMAT (' ',1P, E10.2,' ',E10.2, & & ' tauMax = upper limit(s) on subduction zone shear', & & ' coupling, integrated down-dip (N/m). One value:', & & ' universal. Two values: sea, land.') IF ((tauMax(1) < 0.0D0).OR.(tauMax(2) < 0.0D0)) THEN WRITE (*, 107) WRITE (iUnit8, 107) 107 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' tauMax may not be negative.') CALL Pause() STOP END IF READ (iUnit7, * ) rhoH2O WRITE (iUnit8, 110) rhoH2O 110 FORMAT (' ',11X,1P,E10.3,' rhoH2O = density of groundwater,', & & ' lakes, & oceans') IF (rhoH2O < 0.0D0) THEN WRITE (*, 111) WRITE (iUnit8, 111) 111 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' rhoH2O may not be negative.') CALL Pause() STOP END IF CALL ReadN (iUnit7, iUnit8, 2, & ! input & rhoBar) ! output WRITE (iUnit8, 120) rhoBar(1), rhoBar(2) 120 FORMAT (' ',1P,E10.3,' ',E10.3,' rhoBar = mean density,', & & ' corrected to 0 degrees Kelvin. (crust/mantle)') IF ((rhoBar(1) <= 0.0D0).OR.(rhoBar(2) <= 0.0D0)) THEN WRITE (*, 121) WRITE (iUnit8, 121) 121 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' rhoBar must be positive in each layer.') CALL Pause() STOP END IF READ (iUnit7, * ) rhoAst WRITE (iUnit8, 130) rhoAst 130 FORMAT (' ',11X,1P,E10.3,' rhoAst = density of asthenosphere') IF (rhoAst <= 0.0D0) THEN WRITE (*, 131) WRITE (iUnit8, 131) 131 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' rhoAst must be positive.') CALL Pause() STOP END IF READ (iUnit7, * ) gMean WRITE (iUnit8, 140) gMean 140 FORMAT (' ',11X,1P,E10.3,' gMean = mean gravitational', & & ' acceleration', & & ' (length/s**2)') IF (gMean <= 0.0D0) THEN WRITE (*, 141) WRITE (iUnit8, 141) 141 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' gMean must be positive.') CALL Pause() STOP END IF READ (iUnit7, * ) oneKm WRITE (iUnit8, 150) oneKm 150 FORMAT (' ',11X,1P,E10.3,' oneKm = number of length units', & & ' needed to make 1 kilometer (e.g., 1000. in SI, 1.D5 in cgs)') IF (oneKm <= 0.0D0) THEN WRITE (*, 151) WRITE (iUnit8, 151) 151 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' oneKm must be positive.') CALL Pause() STOP END IF READ (iUnit7, * ) radius WRITE (iUnit8, 155) radius 155 FORMAT (' ',11X,1P,E10.3,' radius = radius of the planet') IF (radius <= 0.0D0) THEN WRITE (*, 156) WRITE (iUnit8, 156) 156 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' radius must be positive.') CALL Pause() STOP END IF CALL ReadN (iUnit7, iUnit8, 2, & ! input & alphaT) ! output WRITE (iUnit8, 160) alphaT(1), alphaT(2) 160 FORMAT (' ',1P,E10.2,' ',E10.2,' alphaT = volumetric thermal', & & ' expansion', & & ' (1/V)*(dV/dT). (crust/mantle)') IF ((alphaT(1) < 0.0D0).OR.(alphaT(2) < 0.0D0)) THEN WRITE (*, 161) WRITE (iUnit8, 161) 161 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' Negative alphaT in either layer is unphysical.') CALL Pause() STOP END IF CALL ReadN (iUnit7, iUnit8, 2, & ! input & conduc) ! output WRITE (iUnit8, 170) conduc(1), conduc(2) 170 FORMAT (' ',1P,E10.2,' ',E10.2,' conduc = thermal conductivity,', & & ' energy/length/s/deg. (crust/mantle)') IF ((conduc(1) <= 0.0D0).OR.(conduc(2) <= 0.0D0)) THEN WRITE (*, 171) WRITE (iUnit8, 171) 171 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' conduc must be positive in each layer.') CALL Pause() STOP END IF CALL ReadN (iUnit7, iUnit8, 2, & ! input & radio) ! output WRITE (iUnit8, 180) radio(1), radio(2) 180 FORMAT (' ',1P,E10.2,' ',E10.2,' radio = radioactive heat ', & & 'production, energy/volume/s. (crust/mantle)') IF ((radio(1) < 0.0D0).OR.(radio(2) < 0.0D0)) THEN WRITE (*, 181) WRITE (iUnit8, 181) 181 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' Negative radio in either layer is unphysical.') CALL Pause() STOP END IF READ (iUnit7, * ) tSurf WRITE (iUnit8, 185) tSurf 185 FORMAT (' ',11X,F10.0,' tSurf = surface temperature, on', & & ' absolute scale (deg. K)') IF (tSurf <= 0.0D0) THEN WRITE (*, 186) WRITE (iUnit8, 186) 186 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' tSurf must be positive.') CALL Pause() STOP END IF CALL ReadN (iUnit7, iUnit8, 2, & ! input & temLim) ! output WRITE (iUnit8, 190) temLim(1), temLim(2) 190 FORMAT (' ',F10.0,' ',F10.0,' temLim = convecting', & & ' temperature (Tmax), on absolute scale. (crust/mantle)') IF ((temLim(1) <= 0.0D0).OR.(temLim(2) <= 0.0D0)) THEN WRITE (*, 191) WRITE (iUnit8, 191) 191 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' temLim must be positive in each layer.') CALL Pause() STOP END IF WRITE (iUnit8, 199) 199 FORMAT (' [OMIT from iXXX.in]' & & /' [OMIT from iXXX.in]', & & ' ========== Tactical parameters (How to reach ', & & 'the solution?) ==========' & & /' [OMIT from iXXX.in]') READ (iUnit7, * ) maxItr WRITE (iUnit8, 200) maxItr 200 FORMAT (' ',11X,I10,' maxItr = maximum iterations within', & & ' velocity solution') IF (maxItr < 1) THEN WRITE (*, 201) WRITE (iUnit8, 201) 201 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' maxItr must be positive.') CALL Pause() STOP END IF READ (iUnit7, * ) okToQt WRITE (iUnit8, 210) okToQt 210 FORMAT (' ',11X,F10.6,' okToQt = acceptable fractional change', & & ' in velocity (stops iteration early)') IF (okToQt < 0.0D0) THEN WRITE (*, 211) WRITE (iUnit8, 211) 211 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' okToQt may not be negative.') CALL Pause() STOP END IF READ (iUnit7, * ) refStr WRITE (iUnit8, 220) refStr 220 FORMAT (' ',11X,1P,E10.2,' refStr = expected mean value of', & & ' shear stress in plate (used to set stiffness limits)') IF (refStr <= 0.0D0) THEN WRITE (*, 221) WRITE (iUnit8, 221) 221 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' refStr must be positive.') CALL Pause() STOP END IF READ (iUnit7, * ) okDelV WRITE (iUnit8, 230) okDelV 230 FORMAT (' ',11X,1P,E10.2,' okDelV = magnitude of velocity', & & ' errors allowed due to finite stiffness' & & /' [OMIT from iXXX.in] ', & & ' (Such errors may appear in such forms as:' & & /' [OMIT from iXXX.in] ', & & ' 1. fictitious basal slip of plate over asthenosphere' & & /' [OMIT from iXXX.in] ', & & ' 2. erroneous convergence/divergence at vertical faults' & & /' [OMIT from iXXX.in] ', & & ' 3. velocity effect of fictitious viscous compliances' & & /' [OMIT from iXXX.in] ', & & ' HOWEVER, VALUES WHICH ARE TOO SMALL WILL CAUSE ILL-CONDITIONED' & & /' [OMIT from iXXX.in] ', & & ' LINEAR SYSTEMS AND STRESS ERR0RS, ', & & 'AND MAY PREVENT CONVERGENCE!)' & &) IF (okDelV <= 0.0D0) THEN WRITE (*, 231) WRITE (iUnit8, 231) 231 FORMAT (/' *** ERR','OR in parameter input file: ***' & & /' okDelV must be positive.') CALL Pause() STOP END IF READ (iUnit7, * ) everyP WRITE (iUnit8, 240) everyP 240 FORMAT (' ',11X,L10,' everyP = should nodal velocities be', & & ' output in every iteration? (for convergence studies)') WRITE (iUnit8, 999) 999 FORMAT (' --------------------------------------------------', & & '-----------------------------') RETURN END SUBROUTINE ReadPm REAL*8 FUNCTION ATan2F (y, x) ! Corrects for problem of two zero arguments. ! N.B. Order of these two arguments is NOT an error; ! this order corresponds to ("motivator", "restrainer"). REAL*8, INTENT(IN) :: y, x IF ((y /= 0.0D0).OR.(x /= 0.0D0)) THEN ATan2F = ATAN2(y, x) ELSE ATan2F = 0.0D0 END IF END FUNCTION ATan2F SUBROUTINE SNodal (phi, theta, & ! input & fpp) ! output ! Calculates all (vector) nodal functions at all integration points along an ! arc-of-great-circle side of any single finite element. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - REAL*8, INTENT(IN) :: phi, theta ! input REAL*8, INTENT(OUT) :: fpp ! output ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - DOUBLE PRECISION fPhi COMMON / FPhis / fPhi DIMENSION fPhi(4, 7) ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - INTEGER j, m REAL*8 cscs, csp, cssn, dd, fp, & & phaij, ppm, rn, sitaj, snc, snp, snsn, & & x1, x2, xn, xx, xyzn, y1, y2, yn, yy, z1, z2, zn, zz DOUBLE PRECISION pp DIMENSION fpp(2, 2, 2, 7), phi(2), theta(2) x1 = SIN(theta(1)) * COS(phi(1)) y1 = SIN(theta(1)) * SIN(phi(1)) z1 = COS(theta(1)) x2 = SIN(theta(2)) * COS(phi(2)) y2 = SIN(theta(2)) * SIN(phi(2)) z2 = COS(theta(2)) xn = x1 + x2 yn = y1 + y2 zn = z1 + z2 xyzn = SQRT(xn * xn + yn * yn + zn * zn) xn = xn / xyzn yn = yn / xyzn zn = zn / xyzn dd = x1 * xn + y1 * yn + z1 * zn DO 800 m = 1, 7 xx = fPhi(1, m) * x1 + fPhi(2, m) * x2 yy = fPhi(1, m) * y1 + fPhi(2, m) * y2 zz = fPhi(1, m) * z1 + fPhi(2, m) * z2 pp = SQRT(xx * xx + yy * yy + zz * zz) xx = xx / pp yy = yy / pp zz = zz / pp sitaj = ACOS(zz) phaij = ATan2F(yy, xx) rn = xx * xn + yy * yn + zz * zn ppm = rn / dd cscs = COS(sitaj) * COS(phaij) cssn = COS(sitaj) * SIN(phaij) snsn = SIN(sitaj) * SIN(phaij) snc = SIN(sitaj) snp = SIN(phaij) csp = COS(phaij) DO 500 j = 1, 2 fp = fPhi(j, m) * ppm fpp(1, 1, j, m) = ( COS(theta(j)) * COS(phi(j)) * cscs & & + COS(theta(j)) * SIN(phi(j)) * cssn & & + SIN(theta(j)) * snc) * fp fpp(2, 1, j, m) = (-SIN(phi(j)) * cscs + COS(phi(j)) * cssn) * fp fpp(1, 2, j, m) = (-COS(theta(j)) * COS(phi(j)) * snp & & + COS(theta(j)) * SIN(phi(j)) * csp) * fp fpp(2, 2, j, m) = ( SIN(phi(j)) * snp & & + COS(phi(j)) * csp) * fp 500 CONTINUE 800 CONTINUE END SUBROUTINE SNodal SUBROUTINE Square (brief, fDip, iUnitT, & ! input & log_strike_adjustments, & & mxBn, mxEl, mxFEl, mxNode, & & mxStar, nFl, nodeF, nodes, & & numEl, numNod, skipBC, radius, wedge, & & xNode, yNode, & ! modify & area, detJ, & ! output & dXS, dYS, dXSP, dYSP, edgeFS, & & edgeTS, fLen, fPFlt, fPSfer, & & fArg, nCond, nodCon, sita, & & checkN, list) ! work ! Check, correct, and complete the geometry of the finite element grid. IMPLICIT NONE ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - LOGICAL, INTENT(IN) :: brief ! input REAL*8, INTENT(IN) :: fDip ! input INTEGER, INTENT(IN) :: iUnitT ! input LOGICAL, INTENT(IN) :: log_strike_adjustments ! input INTEGER, INTENT(IN) :: mxBn, mxEl, mxFEl, mxNode, mxStar, nFl, nodeF, nodes, & ! input & numEl, numNod ! input LOGICAL, INTENT(IN) :: skipBC ! input REAL*8, INTENT(IN) :: radius, wedge ! input REAL*8, INTENT(INOUT) :: xNode, yNode ! modify REAL*8, INTENT(OUT) :: area, detJ, dXS, dYS, dXSP, dYSP ! output LOGICAL, INTENT(OUT) :: edgeFS, edgeTs ! output REAL*8, INTENT(OUT) :: fLen, fPFlt, fPSfer, fArg ! output INTEGER, INTENT(OUT) :: nCond, nodCon ! output REAL*8, INTENT(OUT) :: sita ! output LOGICAL checkN ! work INTEGER list ! work ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - DOUBLE PRECISION fPhi, fPoint, fGauss COMMON / SFault / fPoint COMMON / FPhis / fPhi COMMON / FGList / fGauss DIMENSION fPhi(4, 7), fPoint(7), fGauss(7) ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - CHARACTER*21 obliqu, tag1, tag2, vertic INTEGER i, j, j1, j2, k, kEle, kFault, l, m, & & n, n1, n2, n4, nazi, nazl, nDone, nGood, nInSum, nj1, nj2, & & nl1, nl2, nl3, nl4, nLeft, node, node1, node4, np1, np4, number, nvpart REAL*8 azi, azimut, azl, cosz, & & daz, dazi, dazl, deld, dELon, deld1, deld2, dip1, dip2, dNLat, dx, dy, & & fAngle, phi, phi1, phi2, r, r2, rmax, short, sinz, & & t2, test, theLat, theLon, theta, theta1, theta2, toler, & & x, xmean, xsum, y, ymean, ysum LOGICAL agreed, allOK, found, switch, vert1, vert2 DIMENSION fAngle(2), phi(2), theta(2) DIMENSION area(mxEl), checkN(mxNode), & & detJ(7, mxEl), & & dXS(2, 2, 3, 7, mxEl), dYS(2, 2, 3, 7, mxEl), & & dXSP(3, 7, mxEl), dYSP(3, 7, mxEl), & & edgeFS(2, mxFEl), edgeTS(3, mxEl), fDip(2, mxFEl), & & fLen(mxFEl), & & fPFlt(2, 2, 2, 7, mxFEl), & & fPSfer(2, 2, 3, 7, mxEl), fArg(2, mxFEl), & & list(mxStar), nodCon(mxBn), & & nodeF(4, mxFEl), nodes(3, mxEl), & & sita(7, mxEl), xNode(mxNode), yNode(mxNode) DATA obliqu / '(DIP SLIP IS ALLOWED)' / DATA vertic / '(STRIKE-SLIP ONLY) ' / ! (1) Check that all nodes are connected to at least one ! continuum (triangular) element or fault element: DO 110 i = 1, numNod checkN(i) = .FALSE. 110 CONTINUE DO 130 i = 1, numEl DO 120 j = 1, 3 checkN(nodes(j, i)) = .TRUE. 120 CONTINUE 130 CONTINUE DO 136 i = 1, nFl DO 134 j = 1, 4 checkN(nodeF(j, i)) = .TRUE. 134 CONTINUE 136 CONTINUE allOK = .TRUE. DO 140 i = 1, numNod allOK = allOK.AND.checkN(i) 140 CONTINUE IF (.NOT.allOK) THEN WRITE(iUnitT, 150) 150 FORMAT(' BAD GRID TOPOLOGY: FOLLOWING REAL NODES DO NOT'/ & & ' BELONG TO ANY TRIANGULAR CONTINUUM ELEMENT'/ & & ' OR FAULT ELEMENT:') DO 160 i = 1, numNod IF (.NOT.checkN(i)) WRITE (iUnitT, 155) i 155 FORMAT (' ',43X,I6) 160 CONTINUE CALL Pause() STOP END IF ! (2) Average together the coordinates of all nodes at one "point": DO 410 i = 1, numNod checkN(i) = .FALSE. ! (Means "not yet involved in averaging") 410 CONTINUE DO 490 i = 1, nFl DO 480 j1 = 1, 2 nj1 = nodeF(j1, i) ! (Fault ends are the only places that can have problems.) IF (.NOT.checkN(nj1)) THEN list(1) = nj1 checkN(nj1) = .TRUE. ! Begin list of neighbors with paired node: j2 = 5 - j1 nj2 = nodeF(j2, i) list(2) = nj2 checkN(nj2) = .TRUE. nInSum = 2 ! Find shortest fault connected to either one: dx = xNode(nj1) - xNode(nj2) dy = yNode(nj1) - yNode(nj2) IF (dy > 3.14159265358979D0) dy = dy - 6.28318530717959D0 IF (dy < -3.14159265358979D0) dy = dy + 6.28318530717959D0 dy = dy * SIN(xNode(nj1)) short = SQRT(dx**2 + dy**2) DO 470 k = 1, nFl nl1 = nodeF(1, k) nl2 = nodeF(2, k) nl3 = nodeF(3, k) nl4 = nodeF(4, k) IF ((nj1 == nl1).OR.(nj2 == nl1).OR. & & (nj1 == nl2).OR.(nj2 == nl2).OR. & & (nj1 == nl3).OR.(nj2 == nl3).OR. & & (nj1 == nl4).OR.(nj2 == nl4)) THEN dx = xNode(nl1) - xNode(nl2) dy = yNode(nl1) - yNode(nl2) IF (dy > 3.14159265358979D0) dy = dy - 6.28318530717959D0 IF (dy < -3.14159265358979D0) dy = dy + 6.28318530717959D0 dy = dy * SIN(xNode(nl1)) test = SQRT(dx**2 + dy**2) short = MIN(short, test) END IF 470 CONTINUE ! Collect all corner nodes within 10% of this: toler = short / 10.0D0 t2 = toler**2 DO 471 k = 1, numNod IF (.NOT.checkN(k)) THEN dx = xNode(nj1) - xNode(k) dy = yNode(nj1) - yNode(k) IF (dy > 3.14159265358979D0) dy = dy - 6.28318530717959D0 IF (dy < -3.14159265358979D0) dy = dy + 6.28318530717959D0 dy = dy * SIN(xNode(nj1)) r2 = dx**2 + dy**2 IF (r2 < t2) THEN nInSum = nInSum + 1 IF (nInSum > mxStar) THEN WRITE(iUnitT, 421) 421 FORMAT(/' INCREASE VALUE' & & ,' OF PARAMETER mxStar.') CALL Pause() STOP END IF list(nInSum) = k checkN(k) = .TRUE. END IF END IF 471 CONTINUE ! (Quick EXIT if all nodes in same place) agreed = .TRUE. DO 472 k = 2, nInSum agreed = agreed.AND. & & (xNode(list(k)) == xNode(list(1))).AND. & & (yNode(list(k)) == yNode(list(1))) 472 CONTINUE IF (agreed) GO TO 480 xsum = 0.0D0 ysum = 0.0D0 DO 473 k = 1, nInSum xsum = xsum + xNode(list(k)) ysum = ysum + yNode(list(k)) 473 CONTINUE xmean = xsum / nInSum ymean = ysum / nInSum rmax = 0.0D0 DO 474 k = 1, nInSum r = SQRT((xNode(list(k)) - xmean)**2 + & & (yNode(list(k)) - ymean)**2) rmax = MAX(rmax, r) 474 CONTINUE DO 475 k = 1, nInSum xNode(list(k)) = xmean yNode(list(k)) = ymean 475 CONTINUE IF (.NOT.brief) THEN IF (rmax > 0.0D0) THEN WRITE(iUnitT, 476) nInSum, & & (list(n), n = 1, nInSum) 476 FORMAT(/ & & ' AVERAGING TOGETHER THE POSITIONS OF', & & ' THESE ',I6,' NODES:',(/' ',12I6)) WRITE (iUnitT, 477) rmax 477 FORMAT (' MAXIMUM CORRECTION TO ', & & 'ANY POSITION IS',1P,E10.2/ & & ' YOU ARE RESPONSIBLE FOR ', & & ' DECIDING WHETHER THIS IS A', & & ' SERIOUS ERR0R!') END IF END IF END IF 480 CONTINUE 490 CONTINUE ! (3) Compute derivitives of nodal ! functions at integration points; ! then check for negative areas: CALL Deriv (iUnitT, mxEl, mxNode, nodes, numEl, & ! input & radius, xNode, yNode, & & area, detJ, & ! output & dXS, dYS, dXSP, dYSP, fPSfer, sita) allOK = .TRUE. DO 620 i = 1, numEl DO 610 m = 1, 7 test = area(i) * detJ(m, i) IF (test <= 0.0D0) THEN WRITE(iUnitT, 605) m, i 605 FORMAT(/' EXCESSIVELY DISTORTED ELEMENT LEADS TO ' & & ,'NEGATIVE AREA AT POINT ',I1,' IN ELEMENT ', & & I5) WRITE(iUnitT, 606) area(i), detJ(m, i) 606 FORMAT('AREA = ',1P,E12.4,' detJ: ',0P,F12.6) allOK = .FALSE. END IF 610 CONTINUE 620 CONTINUE IF (.NOT.allOK) THEN CALL Pause() STOP END IF ! (4) Compute lengths of fault elements: DO 750 i = 1, nFl n1 = nodeF(1, i) n2 = nodeF(2, i) theta1 = xNode(n1) theta2 = xNode(n2) phi1 = yNode(n1) phi2 = yNode(n2) fLen(i) = FltLen (phi1, phi2, radius, theta1, theta2) 750 CONTINUE ! (5) Make a list of nodes that are on the boundary and require ! boundary conditions (nodCon); these are in counterclockwise ! order. Also make lists of element sides which contain these ! nodes: edgeTS and edgeFS. nCond = 0 DO 801 i = 1, numNod checkN(i) = .FALSE. 801 CONTINUE DO 802 i = 1, nFl edgeFS(1, i) = .FALSE. edgeFS(2, i) = .FALSE. 802 CONTINUE DO 810 i = 1, numEl DO 809 j = 1, 3 CALL Next (i, j, mxEl, mxFEl, nFl, nodeF, nodes, numEl, & ! input & kFault, kEle) ! output IF (kEle > 0) THEN ! (ordinary interior side) edgeTS(j, i) = .FALSE. ELSE IF (kFault == 0) THEN ! (exterior side) edgeTS(j, i) = .TRUE. n1 = nodes(MOD(j, 3) + 1, i) n2 = nodes(MOD(j + 1, 3) + 1, i) IF (.NOT.checkN(n1)) THEN nCond = nCond + 1 checkN(n1) = .TRUE. END IF IF (.NOT.checkN(n2)) THEN nCond = nCond + 1 checkN(n2) = .TRUE. END IF ELSE ! (triangular element has an exterior fault element ! adjacent to it) edgeTS(j, i) = .FALSE. n1 = nodes(MOD(j, 3) + 1, i) IF (nodeF(2, kFault) == n1) THEN edgeFS(2, kFault) = .TRUE. DO 806 k = 3, 4 n = nodeF(k, kFault) IF (.NOT.checkN(n)) THEN nCond = nCond + 1 checkN(n) = .TRUE. END IF 806 CONTINUE ELSE edgeFS(1, kFault) = .TRUE. DO 808 k = 1, 2 n = nodeF(k, kFault) IF (.NOT.checkN(n)) THEN nCond = nCond + 1 checkN(n) = .TRUE. END IF 808 CONTINUE END IF END IF 809 CONTINUE 810 CONTINUE IF (nCond > mxBn) THEN WRITE(iUnitT, 820) nCond 820 FORMAT(/' Increase array-size mxBn to at least ', I6, & & /' (by adjusting formula) and recompile.') CALL Pause() STOP END IF ! Stop work if no boundary nodes found (global grid): IF (nCond == 0) GO TO 899 ! Begin circuit with lowest-numbered boundary node DO 830 i = 1, numNod IF (checkN(i)) GO TO 831 830 CONTINUE 831 nodCon(1) = i nDone = 1 nLeft = nCond ! Beginning of indefinate loop which traces around the perimeter. ! Each time, it progresses by one of 3 steps: ! -1 node at a time along a triangle side, OR ! -1 node at a time along a fault element side, or ! -by finding another node which shares the same location. ! Beginning of main indefinate loop: 840 node = nodCon(ndone) ! Important: Check that we are not revisiting a node! ! This would mean that there are too many boundary nodes ! to fit in the simply-connected loop, and that there ! are excess boundary nodes somewhere, unconnected! IF (.NOT.checkN(node)) THEN nGood = nDone - 2 WRITE (iUnitT, 841) nGood, nCond 841 FORMAT(/' ERR','OR IN GRID, reported by -Square-:' & & /' BOUNDARY IS NOT SIMPLY-CONNECTED.' & & /' Closed loop of ',I6,' nodes does not' & & /' include all ',I6,' boundary nodes.' & & /' Run command PerimeterTest in -OrbWin-' & & /' for a map of the bad nodes.'/) CALL Pause() STOP END IF IF (nDone > 1) checkN(node) = .FALSE. x = xNode(node) y = yNode(node) ! Look for a triangular element with an external ! side that begins with this node: DO 844 i = 1, numEl DO 842 j = 1, 3 IF (edgeTS(j, i)) THEN n1 = nodes(MOD(j, 3) + 1, i) IF (n1 == node) GO TO 846 END IF 842 CONTINUE 844 CONTINUE GO TO 850 846 n2 = nodes(MOD(j + 1, 3) + 1, i) ! Success by element method: n2 is next boundary node nDone = nDone + 1 IF (nDone <= nCond) nodCon(nDone) = n2 nLeft = nLeft - 1 IF (nLeft > 0) THEN GO TO 840 ELSE GO TO 870 END IF ! Else, look for an adjacent fault element using this node: 850 DO 854 i = 1, nFl IF (edgeFS(1, i)) THEN IF (nodeF(1, i) == node) THEN n2 = nodeF(2, i) GO TO 856 END IF ELSE IF (edgeFS(2, i)) THEN IF (nodeF(3, i) == node) THEN n2 = nodeF(4, i) GO TO 856 END IF END IF 854 CONTINUE GO TO 860 856 nDone = nDone + 1 ! Success by fault method: n2 is next boundary node: IF (nDone <= nCond) nodCon(nDone) = n2 nLeft = nLeft - 1 IF (nLeft > 0) THEN GO TO 840 ELSE GO TO 870 END IF ! Else, look for another exterior corner node at same location: 860 DO 865 i = 1, numNod IF ((i /= node).AND.checkN(i)) THEN IF ( (ABS(xNode(i) - x) < 1.0D-6) .AND. & & (ABS(yNode(i) - y) < 1.0D-6) ) GO TO 867 END IF 865 CONTINUE WRITE(iUnitT, 866) node 866 FORMAT(' BAD GRID TOPOLOGY: WHILE TRACING PERIMETER,'/ & & ' COULD NOT FIND ANY WAY TO CONTINUE FROM NODE ',I6/ & & ' EITHER THROUGH SHARED BOUNDARY ELEMENTS, OR'/ & & ' THROUGH OTHER BOUNDARY NODES SHARING THE SAME ', & & 'POSITION.') CALL Pause() STOP 867 nDone = nDone + 1 ! Success by location method: I is the next boundary node IF (nDone <= nCond) nodCon(nDone) = i nLeft = nLeft - 1 IF (nLeft > 0) GO TO 840 ! End of indefinate loop which traces around perimeter. 870 IF (.NOT.skipBC) THEN WRITE(iUnitT, 880) 880 FORMAT(/ /' Here follows a list, in consecutive order,'/ & & ' of the nodes which define the perimeter'/ & & ' of the model; these nodes require boundary', & & ' conditions:'/' BC# Node ', & & ' Latitude Longitude') DO 890 i = 1, nCond n = nodCon(i) theLon = yNode(n) * 57.2957795130823D0 theLat = 90.0D0 - xNode(n) * 57.2957795130823D0 WRITE(iUnitT, 882) i, n, theLat, theLon 882 FORMAT(' ',2I6,10X,2F10.3) 890 CONTINUE n = nodCon(1) WRITE (iUnitT, 892) n 892 FORMAT(' (Note: Initial node ',I6,' completes the loop,', & & ' but is not listed again.)') END IF 899 CONTINUE ! (6) Survey fault elements and issue warning if any element is of ! mixed type (part strike-slip, and part shallow-dipping): DO 920 i = 1, nFl deld1 = fDip(1, i) - 1.57079632679490D0 deld2 = fDip(2, i) - 1.57079632679490D0 vert1 = ABS(deld1) <= wedge vert2 = ABS(deld2) <= wedge nvpart = 0 IF (vert1) THEN nvpart = nvpart + 1 tag1 = vertic ELSE tag1 = obliqu END IF IF (vert2) THEN nvpart = nvpart + 1 tag2 = vertic ELSE tag2 = obliqu END IF switch = ((nvpart > 0).AND.(nvpart < 2)) IF (switch) THEN dip1 = fDip(1, i) * 57.2957795130823D0 IF (dip1 > 90.0D0) dip1 = dip1 - 180.0D0 dip2 = fDip(2, i) * 57.2957795130823D0 IF (dip2 > 90.0D0) dip2 = dip2 - 180.0D0 WRITE (iUnitT, 905) i, dip1, tag1, dip2, tag2 905 FORMAT(/ /' CAUTION:'/ & & ' FAULT ELEMENT ',I6,' HAS DIPS OF '/ & & ' ',F7.2,' DEGREES ',A21/ & & ' ',F7.2,' DEGREES ',A21/ & & ' WHICH MAKES IT MIXED-MODE.'/ & & ' SUCH ELEMENTS ARE INACCURATE AND NOT RECOMMENDED.'/ & & ' PREFERABLY EACH ELEMENT SHOULD BE OF A SINGLE TYPE.'/ & & ' (REMEMBER, DIP NEED NOT BE CONTINUOUS FROM ONE', & & ' FAULT ELEMENT TO THE NEXT.)') ELSE nvpart = 0 DO 910 m = 1, 7 deld = deld1 * fPhi(1, m) + deld2 * fPhi(2, m) IF (ABS(deld) <= wedge) nvpart = nvpart + 1 910 CONTINUE IF ((nvpart > 0).AND.(nvpart < 7)) THEN IF (nvpart >= 4) THEN WRITE (iUnitT, 912) i, dip1, dip2 912 FORMAT(/ /' CAUTION:'/ & & ' FAULT ELEMENT ',I6,' HAS DIPS OF '/ & & ' ',F7.2,' DEGREES, AND'/ & & ' ',F7.2,' DEGREES'/ & & ' WHICH APPEAR TO MAKE IT STRIKE-SLIP.'/ & & ' HOWEVER, THESE VALUES ARE SUCH THAT DIP-SLIP'/ & & ' IS PERMITTED AT ONE OR MORE INTEGRATION POINTS.'/ & & ' SUCH ELEMENTS ARE INACCURATE AND NOT RECOMMENDED.'/ & & ' PREFERABLY EACH ELEMENT SHOULD BE OF A SINGLE TYPE.'/ & & ' (REMEMBER, DIP NEED NOT BE CONTINUOUS FROM ONE', & & ' FAULT ELEMENT TO THE NEXT.)') ELSE WRITE (iUnitT, 914) i, dip1, dip2 914 FORMAT(/ /' CAUTION:'/ & & ' FAULT ELEMENT ',I6,' HAS DIPS OF '/ & & ' ',F7.2,' DEGREES, AND'/ & & ' ',F7.2,' DEGREES'/ & & ' WHICH APPEAR TO MAKE IT FREE-SLIPPING.'/ & & ' HOWEVER, THESE VALUES ARE SUCH THAT DIP-SLIP'/ & & ' IS PROHIBITED AT ONE OR MORE INTEGRATION POINTS.'/ & & ' SUCH ELEMENTS ARE INACCURATE AND NOT RECOMMENDED.'/ & & ' PREFERABLY EACH ELEMENT SHOULD BE OF A SINGLE TYPE.'/ & & ' (REMEMBER, DIP NEED NOT BE CONTINUOUS FROM ONE', & & ' FAULT ELEMENT TO THE NEXT.)') END IF END IF END IF 920 CONTINUE ! (7) Calculate fault argument (in radians, measured counterclockwise ! from +Theta = South) at each end of each fault element. DO 1000 i = 1, nFl n1 = nodeF(1, i) n2 = nodeF(2, i) theta(1) = xNode(n1) theta(2) = xNode(n2) phi(1) = yNode(n1) phi(2) = yNode(n2) CALL FAngls(phi, theta, & ! input & fAngle) ! output DO 900 j = 1, 2 fArg(j, i) = fAngle(j) 900 CONTINUE 1000 CONTINUE ! (8) Survey strike-slip (vertical) faults to check for conflicts in ! argument that would lock the fault: IF (log_strike_adjustments) WRITE(iUnitT, 1001) 1001 FORMAT(/ /' The following tightly-connected pairs of strike-slip' & & /' fault elements had their azimuths averaged at the' & & /' connection point for purposes of computing the' & & /' constraint on the direction of strike-slip:' & & / /' Fault#1 Fault#2 Node#A Node#B ', & & ' Latitude Longitude Azim#1 Azim#2 Azimuth' & & /' ----------------------------------------', & & '--------------------------------------------------') ! Loop on all fault elements (I): DO 2000 i = 1, nFl ! Loop on 2 terminal node pairs, 1-4, 2-3 (J = 1 or 2): DO 1900 j = 1, 2 ! Dip must be within "wedge" of vertical for constraint: IF (ABS(fDip(j, i) - 1.57079632679490D0) <= wedge) THEN nazi = j n1 = j IF(j == 1) THEN n4 = 4 ELSE n4 = 3 END IF node1 = nodeF(n1, i) node4 = nodeF(n4, i) ! No constraint applied where a fault ends: IF (node1 /= node4) THEN ! Endpoint pairs must be checkEd for duplication: ! Look for other strike-slip faults sharing this ! pair of nodes, at either end: found = .FALSE. DO 1600 l = 1, nFl IF (l /= i) THEN IF (ABS(fDip(1, l) - 1.57079632679490D0) <= wedge) THEN IF (((node1 == nodeF(1, l)).AND. & & (node4 == nodeF(4, l))).OR. & & ((node1 == nodeF(4, l)).AND. & & (node4 == nodeF(1, l)))) THEN found = .TRUE. number = l nazl = 1 GO TO 1601 END IF END IF IF (ABS(fDip(2, l) - 1.57079632679490D0) <= wedge) THEN IF (((node1 == nodeF(2, l)).AND. & & (node4 == nodeF(3, l))).OR. & & ((node1 == nodeF(3, l)).AND. & & (node4 == nodeF(2, l)))) THEN found = .TRUE. number = l nazl = 2 GO TO 1601 END IF END IF END IF 1600 CONTINUE ! Don't worry if this pair already checkEd! 1601 IF (found.AND.(number > i)) THEN ! Average arguments together (avoid cycle shifts): IF(nazi == nazl) THEN azl = fArg(nazl, number) + 3.14159265358979D0 ELSE azl = fArg(nazl, number) END IF azi = fArg(nazi, i) cosz = 0.5D0 * (COS(azi) + COS(azl)) sinz = 0.5D0 * (SIN(azi) + SIN(azl)) azimut = ATan2F(sinz, cosz) fArg(nazi, i) = azimut IF(nazl == nazi) THEN fArg(nazl, number) = azimut - 3.14159265358979D0 ELSE fArg(nazl, number) = azimut END IF ! Print a warning: dazi = azi * 57.2957795130823D0 dazl = azl * 57.2957795130823D0 daz = azimut * 57.2957795130823D0 np1 = node1 np4 = node4 dELon = 57.2957795130823D0 * yNode(node1) dNLat = 90.0D0 - 57.2957795130823D0 * xNode(node1) IF (log_strike_adjustments) WRITE (iUnitT, 1610) & & i, number, np1, np4, & & dNLat, dELon, dazi, dazl, daz 1610 FORMAT(' ',I7,3X,I7,3X, & & I7,3X,I7,3X, & & 2X,F8.3,1X,F9.3, & & 4X,F6.1,4X,F6.1,4X,F6.1) END IF ! ^End block which looks for constraints END IF ! ^End block which checks for distinct node numbers END IF ! ^End block which checks for dip of over 75 degrees 1900 CONTINUE ! ^End loop on 2 node pairs in fault element 2000 CONTINUE ! (9) Calculate nodal functions at integration points on faults: CALL FNodal (mxFEl, & ! input & mxNode, nFl, nodeF, xNode, yNode, & & fPFlt) ! output IF (.NOT. brief) WRITE (iUnitT, 9999) 9999 FORMAT (' -------------------------------------------------------------------------------') END SUBROUTINE Square SUBROUTINE ThetaPhi_2_plate (iUnitT, & & nPBnd, nBoundaryPoints, & & nPlate, pLat, pLon, & & theta, phi, & ! INTENT(IN) & iPlate) ! INTENT(OUT) ! Assigns an integer plate ID# to one surface point ! with coordinates (theta = colatitude, in radians, ! phi = longitude, in radians). ! NOTE: This subprogram should NOT be called to assign ! F-E nodes to plates if those nodes might belong to ! fault elements along plate edges! Results would be purely ! random, and would not consider topological connectivity ! to an adjacent plate. Instead, use subprogram Assign ! from -Shells- (versions of 2006+). USE DSphere ! from Peter Bird's file DSphere.f90 IMPLICIT NONE INTEGER, INTENT(IN) :: iUnitT INTEGER, INTENT(IN) :: nPBnd INTEGER, DIMENSION(nPlate), INTENT(IN) :: nBoundaryPoints INTEGER, INTENT(IN) :: nPlate REAL*8, DIMENSION(nPlate, nPBnd), INTENT(IN):: pLat, pLon REAL*8, INTENT(IN) :: theta, phi INTEGER, INTENT(OUT) :: iPlate INTEGER :: iP, j, j2, nEnd, nPoint REAL*8 :: a1, a2, a3, aa, ab1, ab2, ab3, angle, ao, b1, b2, b3, bb, bo, dAngle, & & oxyz, sTheta, tangl, xo, xPoint, yo, yPoint, zo ! PB2002 model of Bird [2003; G**3]; ! Already has plate "names" and "omega" vectors in ! main program (DATA statements); ! must also have digitised plate ! outlines in arrays pLat and pLon, ! presumably read from file "PB2002_plates.dig". ! That is, this routine will not read any file. xo = COS(phi) * SIN(theta) yo = SIN(phi) * SIN(theta) zo = COS(theta) oxyz = xo * xo + yo * yo + zo * zo oxyz = SQRT(oxyz) xo = xo / oxyz yo = yo / oxyz zo = zo / oxyz nPoint = 0 angle = 0.0D0 iPlate = 0 DO 500 iP = 1, nPlate tangl = 0.0D0 nEnd = nBoundaryPoints(ip) DO 300 j = 1, nEnd j2 = j + 1 IF(j == nEnd) THEN j2 = 1 END IF a1 = COS(pLon(ip, j)) * COS(pLat(ip, j)) a2 = SIN(pLon(ip, j)) * COS(pLat(ip, j)) a3 = SIN(pLat(ip, j)) b1 = COS(pLon(ip, j2)) * COS(pLat(ip, j2)) b2 = SIN(pLon(ip, j2)) * COS(pLat(ip, j2)) b3 = SIN(pLat(ip, j2)) ao = xo * a1 + yo * a2 + zo * a3 bo = xo * b1 + yo * b2 + zo * b3 a1 = a1 / ao a2 = a2 / ao a3 = a3 / ao b1 = b1 / bo b2 = b2 / bo b3 = b3 / bo a1 = a1 - xo a2 = a2 - yo a3 = a3 - zo b1 = b1 - xo b2 = b2 - yo b3 = b3 - zo aa = SQRT(a1 * a1 + a2 * a2 + a3 * a3) bb = SQRT(b1 * b1 + b2 * b2 + b3 * b3) ab1 = a2 * b3 - a3 * b2 ab2 = a3 * b1 - a1 * b3 ab3 = a1 * b2 - a2 * b1 sTheta = (ab1 * xo + ab2 * yo + ab3 * zo) / (aa * bb) ! prevent stupid abends due to imprecision: sTheta = MAX(-1.0D0, MIN(1.0D0, stheta)) tangl = tangl + ASIN(stheta) 300 CONTINUE dAngle = tangl - Pi IF(dAngle >= 0.0001D0) THEN nPoint = nPoint + 1 iPlate = iP END IF 500 CONTINUE IF(npoint >= 3) THEN xpoint = 90.0D0 - theta * degrees_per_radian ypoint = phi * degrees_per_radian WRITE(iUnitT, 505) xpoint, ypoint 505 FORMAT(' POINT ', 2F10.3, ' WAS FOUND IN MORE THAN TWO PLATES; SOMETHING IS WRONG !!!!') CALL Pause() STOP END IF IF (iPlate == 0) THEN xPoint = 90.0D0 - theta * degrees_per_radian yPoint = phi * degrees_per_radian WRITE(iUnitT, 600) xpoint, ypoint CALL Pause() STOP 600 FORMAT('THE POINT ', 2F10.2, ' DOES NOT BELONG TO ANY PLATE !!!!') CALL Pause() STOP END IF END SUBROUTINE ThetaPhi_2_plate SUBROUTINE Tractor(iUnitQ, iUnitT, nPlates, nPoints, & & slab_q, phi_list, theta_list, whichP, & ! INTENT(IN) & basal_shear_tractions) ! INTENT(OUT) ! Requests file name of an existing torque report ! (including traction pole vectors for each plate, ! created by a previous experiment with -Shells-, ! usually one that had trHMax = 0.0 and extra internal ! velocity boundary conditions for each slabless plate). ! Reads this file, extracts the traction pole vectors, ! and uses them to precompute basal_shear_tractions shear tractions ! at each point specified in (theta_list, phi_list) = ! (colatitude, longitude) in radians. These points must ! already be identified as to plate affinity, in whichP. ! For further clarification of "traction pole vectors" ! see subprogram -TWIST- in SHELLS. ! Derived from subprogram Tract of Shells (summer 2006), ! but modified to free-format, and variable names changed ! to drop previous implication that the list of points ! to be processed is a list of node locations. USE DSphere ! in Peter Bird's file DSphere.f90 IMPLICIT NONE INTEGER, INTENT(IN) :: iUnitQ, iUnitT, nPlates, nPoints LOGICAL, DIMENSION(nPlates), INTENT(IN) :: slab_q REAL*8, DIMENSION(nPoints), INTENT(IN) :: theta_list, phi_list INTEGER, DIMENSION(nPoints), INTENT(IN) :: whichP REAL*8, DIMENSION(2, nPoints), INTENT(OUT):: basal_shear_tractions CHARACTER*132 :: line INTEGER :: i, ios, iPlate, j LOGICAL, DIMENSION(:), ALLOCATABLE :: tpread REAL*8 :: equat, lat, length, lon, t, tequat REAL*8, DIMENSION(3) :: tvec, uPhi, uTheta, uvec REAL*8, DIMENSION(:, :), ALLOCATABLE :: tpvecs ALLOCATE ( tpvecs(3, nPlates) ) ALLOCATE ( tpread(nPlates) ) ! Zero whole array; advisable because some plates may not appear in the report. DO 30 j = 1, nPlates DO 20 i = 1, 3 tpvecs(i, j) = 0.0D0 20 CONTINUE tpread(j) = .FALSE. 30 CONTINUE WRITE(iUnitT, 35) iUnitQ 35 FORMAT(/' ATTEMPTING TO READ TORQUE REPORT ON UNIT', I3/) ! Waste first 6 lines (titles & 2 blanks & header) of torque file. DO 40 i = 1, 6 READ (iUnitQ, "(A)") line 40 CONTINUE ! Loop on plates in report (up to nPlates for whole-Earth model): DO 100 j = 1, nPlates READ(iUnitQ, * , IOSTAT = ios) IF (ios == -1) GO TO 101 READ(iUnitQ, "(8X,I6)", IOSTAT = ios) iPlate IF (ios == -1) GO TO 101 ! Waste 23 more lines of each plate report DO 50 i = 1, 23 READ(iUnitQ, * , IOSTAT = ios) IF (ios == -1) GO TO 101 50 CONTINUE READ(iUnitQ, "(56X,ES10.3,2F10.2)") t, lon, lat ! T is magnitude, in Pa, at location 90 deg. from (LON, LAT). tpvecs(1, iPlate) = t * COS(lat * radians_per_degree) * COS(lon * radians_per_degree) tpvecs(2, iPlate) = t * COS(lat * radians_per_degree) * SIN(lon * radians_per_degree) tpvecs(3, iPlate) = t * SIN(lat * radians_per_degree) tpread(iPlate) = .TRUE. ! Waste 14 lines to get past the "=======" at the bottom of ! each torque report: DO 60 i = 1, 14 READ(iUnitQ, * , IOSTAT = ios) IF (ios == -1) GO TO 101 60 CONTINUE 100 CONTINUE 101 CLOSE(iUnitQ) DO 200 i = 1, nPoints iPlate = whichP(i) IF (slab_q(iPlate)) THEN ! no need for inferred basal_shear_tractions-strength traction: basal_shear_tractions(1, i) = 0.0D0 basal_shear_tractions(2, i) = 0.0D0 ELSE IF (tpread(iPlate)) THEN ! Uvec is unit vector to node location: uvec(1) = SIN(theta_list(i)) * COS(phi_list(i)) uvec(2) = SIN(theta_list(i)) * SIN(phi_list(i)) uvec(3) = COS(theta_list(i)) ! Tvec is cross-product with traction pole vector: tvec(1) = tpvecs(2, iPlate) * uvec(3) - & & tpvecs(3, iPlate) * uvec(2) tvec(2) = tpvecs(3, iPlate) * uvec(1) - & & tpvecs(1, iPlate) * uvec(3) tvec(3) = tpvecs(1, iPlate) * uvec(2) - & & tpvecs(2, iPlate) * uvec(1) t = SQRT(tvec(1)**2 + tvec(2)**2 + tvec(3)**2) ! Unit vectors at this site (NOT a pole): uPhi(1) = -uvec(2) uPhi(2) = uvec(1) equat = SIN(theta_list(i)) uPhi(1) = uPhi(1) / equat uPhi(2) = uPhi(2) / equat uPhi(3) = 0.0D0 tequat = uvec(3) uTheta(3) = -equat uTheta(1) = tequat * uvec(1) / equat uTheta(2) = tequat * uvec(2) / equat length = SQRT(uTheta(1)**2 + uTheta(2)**2 + uTheta(3)**2) IF (length > 0.0D0) THEN uTheta(1) = uTheta(1) / length uTheta(2) = uTheta(2) / length uTheta(3) = uTheta(3) / length END IF ! Horizontal components of shear traction: basal_shear_tractions(1, i) = tvec(1) * uTheta(1) + tvec(2) * uTheta(2) + & & tvec(3) * uTheta(3) basal_shear_tractions(2, i) = tvec(1) * uPhi(1) + tvec(2) * uPhi(2) + & & tvec(3) * uPhi(3) ELSE basal_shear_tractions(1, i) = 0.0D0 basal_shear_tractions(2, i) = 0.0D0 END IF END IF 200 CONTINUE DEALLOCATE ( tpread ) DEALLOCATE ( tpvecs ) END SUBROUTINE Tractor SUBROUTINE Unit (anyVec) ! Converts any 3-component vector to a unit vector. REAL*8, DIMENSION(3), INTENT(INOUT) :: anyVec REAL*8 :: r2, size r2 = anyVec(1) * anyVec(1) + anyVec(2) * anyVec(2) + anyVec(3) * anyVec(3) IF (r2 > 0.0D0) THEN size = 1.0D0 / SQRT(r2) anyVec(1) = anyVec(1) * size anyVec(2) = anyVec(2) * size anyVec(3) = anyVec(3) * size ELSE anyVec(1) = 1.0D0 anyVec(2) = 0.0D0 anyVec(3) = 0.0D0 END IF END SUBROUTINE Unit END PROGRAM OrbScore2